modbc.F90 27 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649
  1. !
  2. ! $Id: modbc.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_Boundary.
  25. !>
  26. !> Contains subroutines to calculate the boundary conditions on the child grids from their
  27. !> parent grids.
  28. !
  29. module Agrif_Boundary
  30. !
  31. use Agrif_Interpolation
  32. !
  33. implicit none
  34. !
  35. contains
  36. !
  37. !===================================================================================================
  38. ! subroutine Agrif_CorrectVariable
  39. !
  40. !> subroutine to calculate the boundary conditions on a fine grid
  41. !---------------------------------------------------------------------------------------------------
  42. subroutine Agrif_CorrectVariable ( parent, child, pweight, weight, procname )
  43. !---------------------------------------------------------------------------------------------------
  44. type(Agrif_Variable), pointer :: parent !< Variable on the parent grid
  45. type(Agrif_Variable), pointer :: child !< Variable on the child grid
  46. logical :: pweight !< Indicates if weight is used for the time interpolation
  47. real :: weight !< Coefficient for the time interpolation
  48. procedure() :: procname !< Data recovery procedure
  49. !
  50. type(Agrif_Grid) , pointer :: Agrif_Child_Gr, Agrif_Parent_Gr
  51. type(Agrif_Variable), pointer :: root_var ! Variable on the root grid
  52. integer :: nbdim ! Number of dimensions of the grid variable
  53. integer :: n
  54. integer, dimension(6) :: lb_child ! Index of the first point inside the domain for
  55. ! the child grid variable
  56. integer, dimension(6) :: lb_parent ! Index of the first point inside the domain for
  57. ! the parent grid variable
  58. integer, dimension(6) :: ub_child ! Upper bound on the child grid
  59. integer, dimension(6) :: nb_child ! Number of cells for child
  60. integer, dimension(6) :: posvartab_child ! Position of the variable on the cell
  61. integer, dimension(6) :: loctab_child ! Indicates if the child grid has a common border
  62. ! with the root grid
  63. real, dimension(6) :: s_child, s_parent ! Positions of the parent and child grids
  64. real, dimension(6) :: ds_child, ds_parent ! Space steps of the parent and child grids
  65. !
  66. call PreProcessToInterpOrUpdate( parent, child, &
  67. nb_child, ub_child, &
  68. lb_child, lb_parent, &
  69. s_child, s_parent, &
  70. ds_child, ds_parent, nbdim, interp=.true.)
  71. root_var => child % root_var
  72. Agrif_Child_Gr => Agrif_Curgrid
  73. Agrif_Parent_Gr => Agrif_Curgrid % parent
  74. !
  75. loctab_child(:) = 0
  76. posvartab_child(1:nbdim) = root_var % posvar(1:nbdim)
  77. !
  78. do n = 1,nbdim
  79. !
  80. select case(root_var % interptab(n))
  81. !
  82. case('x') ! x DIMENSION
  83. !
  84. if (Agrif_Curgrid % NearRootBorder(1)) loctab_child(n) = -1
  85. if (Agrif_Curgrid % DistantRootBorder(1)) loctab_child(n) = -2
  86. if ((Agrif_Curgrid % NearRootBorder(1)) .AND. &
  87. (Agrif_Curgrid % DistantRootBorder(1))) loctab_child(n) = -3
  88. !
  89. case('y') ! y DIMENSION
  90. !
  91. if (Agrif_Curgrid % NearRootBorder(2)) loctab_child(n) = -1
  92. if (Agrif_Curgrid % DistantRootBorder(2)) loctab_child(n) = -2
  93. if ((Agrif_Curgrid % NearRootBorder(2)) .AND. &
  94. (Agrif_Curgrid % DistantRootBorder(2))) loctab_child(n) = -3
  95. !
  96. case('z') ! z DIMENSION
  97. !
  98. if (Agrif_Curgrid % NearRootBorder(3)) loctab_child(n) = -1
  99. if (Agrif_Curgrid % DistantRootBorder(3)) loctab_child(n) = -2
  100. if ((Agrif_Curgrid % NearRootBorder(3)) .AND. &
  101. (Agrif_Curgrid % DistantRootBorder(3))) loctab_child(n) = -3
  102. !
  103. case('N') ! No space DIMENSION
  104. !
  105. posvartab_child(n) = 1
  106. loctab_child(n) = -3
  107. !
  108. end select
  109. !
  110. enddo
  111. !
  112. call Agrif_Correctnd(parent, child, pweight, weight, &
  113. lb_child(1:nbdim), lb_parent(1:nbdim), &
  114. nb_child(1:nbdim), posvartab_child(1:nbdim), &
  115. loctab_child(1:nbdim), &
  116. s_child(1:nbdim), s_parent(1:nbdim), &
  117. ds_child(1:nbdim),ds_parent(1:nbdim), nbdim, procname )
  118. !---------------------------------------------------------------------------------------------------
  119. end subroutine Agrif_CorrectVariable
  120. !===================================================================================================
  121. !
  122. !===================================================================================================
  123. ! subroutine Agrif_Correctnd
  124. !
  125. !> calculates the boundary conditions for a nD grid variable on a fine grid by using
  126. !> a space and time interpolations; it is called by the #Agrif_CorrectVariable procedure
  127. !---------------------------------------------------------------------------------------------------
  128. subroutine Agrif_Correctnd ( parent, child, pweight, weight, &
  129. pttab_child, pttab_Parent, &
  130. nbtab_Child, posvartab_Child, loctab_Child, &
  131. s_Child, s_Parent, ds_Child, ds_Parent, &
  132. nbdim, procname )
  133. !---------------------------------------------------------------------------------------------------
  134. #if defined AGRIF_MPI
  135. include 'mpif.h'
  136. #endif
  137. !
  138. TYPE(Agrif_Variable), pointer :: parent !< Variable on the parent grid
  139. TYPE(Agrif_Variable), pointer :: child !< Variable on the child grid
  140. LOGICAL :: pweight !< Indicates if weight is used for the temporal interpolation
  141. REAL :: weight !< Coefficient for the temporal interpolation
  142. INTEGER, DIMENSION(nbdim) :: pttab_child !< Index of the first point inside the domain for the parent grid variable
  143. INTEGER, DIMENSION(nbdim) :: pttab_Parent !< Index of the first point inside the domain for the child grid variable
  144. INTEGER, DIMENSION(nbdim) :: nbtab_Child !< Number of cells of the child grid
  145. INTEGER, DIMENSION(nbdim) :: posvartab_Child !< Position of the grid variable (1 or 2)
  146. INTEGER, DIMENSION(nbdim) :: loctab_Child !< Indicates if the child grid has a common border with the root grid
  147. REAL , DIMENSION(nbdim) :: s_Child, s_Parent !< Positions of the parent and child grids
  148. REAL , DIMENSION(nbdim) :: ds_Child, ds_Parent !< Space steps of the parent and child grids
  149. INTEGER :: nbdim !< Number of dimensions of the grid variable
  150. procedure() :: procname !< Data recovery procedure
  151. !
  152. INTEGER,DIMENSION(6) :: type_interp ! Type of interpolation (linear, spline,...)
  153. INTEGER,DIMENSION(6,6) :: type_interp_bc ! Type of interpolation (linear, spline,...)
  154. INTEGER,DIMENSION(nbdim,2,2) :: childarray
  155. INTEGER,DIMENSION(nbdim,2) :: lubglob
  156. INTEGER :: kindex ! Index used for safeguard and time interpolation
  157. INTEGER,DIMENSION(nbdim,2,2) :: indtab ! Arrays indicating the limits of the child
  158. INTEGER,DIMENSION(nbdim,2,2) :: indtruetab ! grid variable where boundary conditions are
  159. INTEGER,DIMENSION(nbdim,2,2,nbdim) :: ptres,ptres2 ! calculated
  160. INTEGER,DIMENSION(nbdim) :: coords
  161. INTEGER :: i, nb, ndir
  162. INTEGER :: n, sizetab
  163. INTEGER :: ibeg, iend
  164. INTEGER :: i1,i2,j1,j2,k1,k2,l1,l2,m1,m2,n1,n2
  165. REAL :: c1t,c2t ! Coefficients for the time interpolation (c2t=1-c1t)
  166. #if defined AGRIF_MPI
  167. !
  168. INTEGER, DIMENSION(nbdim) :: lower, upper
  169. INTEGER, DIMENSION(nbdim) :: ltab, utab
  170. !
  171. #endif
  172. !
  173. type_interp_bc = child % root_var % type_interp_bc
  174. coords = child % root_var % coords
  175. !
  176. ibeg = child % bcinf
  177. iend = child % bcsup
  178. !
  179. indtab(1:nbdim,2,1) = pttab_child(1:nbdim) + nbtab_child(1:nbdim) + ibeg
  180. indtab(1:nbdim,2,2) = indtab(1:nbdim,2,1) + ( iend - ibeg )
  181. indtab(1:nbdim,1,1) = pttab_child(1:nbdim) - iend
  182. indtab(1:nbdim,1,2) = pttab_child(1:nbdim) - ibeg
  183. WHERE (posvartab_child(1:nbdim) == 2)
  184. indtab(1:nbdim,1,1) = indtab(1:nbdim,1,1) - 1
  185. indtab(1:nbdim,1,2) = indtab(1:nbdim,1,2) - 1
  186. END WHERE
  187. !
  188. call Agrif_get_var_global_bounds(child,lubglob,nbdim)
  189. !
  190. indtruetab(1:nbdim,1,1) = max(indtab(1:nbdim,1,1), lubglob(1:nbdim,1))
  191. indtruetab(1:nbdim,1,2) = max(indtab(1:nbdim,1,2), lubglob(1:nbdim,1))
  192. indtruetab(1:nbdim,2,1) = min(indtab(1:nbdim,2,1), lubglob(1:nbdim,2))
  193. indtruetab(1:nbdim,2,2) = min(indtab(1:nbdim,2,2), lubglob(1:nbdim,2))
  194. !
  195. do nb = 1,nbdim
  196. do ndir = 1,2
  197. !
  198. if (loctab_child(nb) /= (-ndir) .AND. loctab_child(nb) /= -3) then
  199. !
  200. do n = 1,2
  201. ptres(nb,n,ndir,nb) = indtruetab(nb,ndir,n)
  202. enddo
  203. !
  204. do i = 1,nbdim
  205. !
  206. if (i /= nb) then
  207. !
  208. if (loctab_child(i) == -1 .OR. loctab_child(i) == -3) then
  209. ptres(i,1,ndir,nb) = pttab_child(i)
  210. else
  211. ptres(i,1,ndir,nb) = indtruetab(i,1,1)
  212. endif
  213. if (loctab_child(i) == -2 .OR. loctab_child(i) == -3) then
  214. if (posvartab_child(i) == 1) then
  215. ptres(i,2,ndir,nb) = pttab_child(i) + nbtab_child(i)
  216. else
  217. ptres(i,2,ndir,nb) = pttab_child(i) + nbtab_child(i) - 1
  218. endif
  219. else
  220. ptres(i,2,ndir,nb) = indtruetab(i,2,2)
  221. endif
  222. !
  223. endif
  224. !
  225. enddo
  226. !
  227. #if defined AGRIF_MPI
  228. call Agrif_get_var_bounds_array(child,lower,upper,nbdim)
  229. do i = 1,nbdim
  230. !
  231. Call Agrif_GetLocalBoundaries(ptres(i,1,ndir,nb), ptres(i,2,ndir,nb), &
  232. coords(i), lower(i), upper(i), ltab(i), utab(i) )
  233. ptres2(i,1,ndir,nb) = max(ltab(i),lower(i))
  234. ptres2(i,2,ndir,nb) = min(utab(i),upper(i))
  235. if ((i == nb) .AND. (ndir == 1)) then
  236. ptres2(i,2,ndir,nb) = max(utab(i),lower(i))
  237. elseif ((i == nb) .AND. (ndir == 2)) then
  238. ptres2(i,1,ndir,nb) = min(ltab(i),upper(i))
  239. endif
  240. !
  241. enddo
  242. #else
  243. ptres2(:,:,ndir,nb) = ptres(:,:,ndir,nb)
  244. #endif
  245. endif
  246. !
  247. enddo ! ndir = 1,2
  248. enddo ! nb = 1,nbdim
  249. !
  250. if ( child % interpIndex /= Agrif_Curgrid % parent % ngridstep .OR. &
  251. child % Interpolationshouldbemade ) then
  252. !
  253. ! Space interpolation
  254. !
  255. kindex = 1
  256. !
  257. do nb = 1,nbdim
  258. type_interp = type_interp_bc(nb,:)
  259. do ndir = 1,2
  260. !
  261. if (loctab_child(nb) /= (-ndir) .AND. loctab_child(nb) /= -3) then
  262. !
  263. call Agrif_InterpnD(type_interp, parent, child, &
  264. ptres(1:nbdim,1,ndir,nb), &
  265. ptres(1:nbdim,2,ndir,nb), &
  266. pttab_child(1:nbdim), &
  267. pttab_Parent(1:nbdim), &
  268. s_Child(1:nbdim), s_Parent(1:nbdim), &
  269. ds_Child(1:nbdim),ds_Parent(1:nbdim), &
  270. NULL(), .FALSE., nbdim, &
  271. childarray, &
  272. child%memberin(nb,ndir), .TRUE., procname, coords(nb), ndir)
  273. child % childarray(1:nbdim,:,:,nb,ndir) = childarray
  274. if (.not. child%interpolationshouldbemade) then
  275. !
  276. ! Safeguard of the values of the grid variable (at times n and n+1 on the parent grid)
  277. !
  278. sizetab = 1
  279. do i = 1,nbdim
  280. sizetab = sizetab * (ptres2(i,2,ndir,nb)-ptres2(i,1,ndir,nb)+1)
  281. enddo
  282. call saveAfterInterp(child,ptres2(:,:,ndir,nb),kindex,sizetab,nbdim)
  283. !
  284. endif
  285. !
  286. endif
  287. !
  288. enddo ! ndir = 1,2
  289. enddo ! nb = 1,nbdim
  290. !
  291. child % interpIndex = Agrif_Curgrid % parent % ngridstep
  292. !
  293. endif
  294. !
  295. if (.not. child%interpolationshouldbemade) then
  296. !
  297. ! Calculation of the coefficients c1t and c2t for the temporary interpolation
  298. !
  299. if (pweight) then
  300. c1t = weight
  301. else
  302. c1t = (REAL(Agrif_Nbstepint()) + 1.) / Agrif_Rhot()
  303. endif
  304. c2t = 1. - c1t
  305. !
  306. ! Time interpolation
  307. !
  308. kindex = 1
  309. !
  310. do nb = 1,nbdim
  311. do ndir = 1,2
  312. if (loctab_child(nb) /= (-ndir) .AND. loctab_child(nb) /= -3) then
  313. Call timeInterpolation(child,ptres2(:,:,ndir,nb),kindex,c1t,c2t,nbdim)
  314. endif
  315. enddo
  316. enddo
  317. !
  318. endif
  319. !
  320. do nb = 1,nbdim
  321. do ndir = 1,2
  322. if ( (loctab_child(nb) /= (-ndir)) .AND. (loctab_child(nb) /= -3) .AND. child%memberin(nb,ndir) ) then
  323. select case(nbdim)
  324. case(1)
  325. i1 = child % childarray(1,1,2,nb,ndir)
  326. i2 = child % childarray(1,2,2,nb,ndir)
  327. call procname(parray1(i1:i2), &
  328. i1,i2, .FALSE.,coords(nb),ndir)
  329. case(2)
  330. i1 = child % childarray(1,1,2,nb,ndir)
  331. i2 = child % childarray(1,2,2,nb,ndir)
  332. j1 = child % childarray(2,1,2,nb,ndir)
  333. j2 = child % childarray(2,2,2,nb,ndir)
  334. call procname(parray2(i1:i2,j1:j2), &
  335. i1,i2,j1,j2, .FALSE.,coords(nb),ndir)
  336. case(3)
  337. i1 = child % childarray(1,1,2,nb,ndir)
  338. i2 = child % childarray(1,2,2,nb,ndir)
  339. j1 = child % childarray(2,1,2,nb,ndir)
  340. j2 = child % childarray(2,2,2,nb,ndir)
  341. k1 = child % childarray(3,1,2,nb,ndir)
  342. k2 = child % childarray(3,2,2,nb,ndir)
  343. call procname(parray3(i1:i2,j1:j2,k1:k2), &
  344. i1,i2,j1,j2,k1,k2, .FALSE.,coords(nb),ndir)
  345. case(4)
  346. i1 = child % childarray(1,1,2,nb,ndir)
  347. i2 = child % childarray(1,2,2,nb,ndir)
  348. j1 = child % childarray(2,1,2,nb,ndir)
  349. j2 = child % childarray(2,2,2,nb,ndir)
  350. k1 = child % childarray(3,1,2,nb,ndir)
  351. k2 = child % childarray(3,2,2,nb,ndir)
  352. l1 = child % childarray(4,1,2,nb,ndir)
  353. l2 = child % childarray(4,2,2,nb,ndir)
  354. call procname(parray4(i1:i2,j1:j2,k1:k2,l1:l2), &
  355. i1,i2,j1,j2,k1,k2,l1,l2, .FALSE.,coords(nb),ndir)
  356. case(5)
  357. i1 = child % childarray(1,1,2,nb,ndir)
  358. i2 = child % childarray(1,2,2,nb,ndir)
  359. j1 = child % childarray(2,1,2,nb,ndir)
  360. j2 = child % childarray(2,2,2,nb,ndir)
  361. k1 = child % childarray(3,1,2,nb,ndir)
  362. k2 = child % childarray(3,2,2,nb,ndir)
  363. l1 = child % childarray(4,1,2,nb,ndir)
  364. l2 = child % childarray(4,2,2,nb,ndir)
  365. m1 = child % childarray(5,1,2,nb,ndir)
  366. m2 = child % childarray(5,2,2,nb,ndir)
  367. call procname(parray5(i1:i2,j1:j2,k1:k2,l1:l2,m1:m2), &
  368. i1,i2,j1,j2,k1,k2,l1,l2,m1,m2, .FALSE.,coords(nb),ndir)
  369. case(6)
  370. i1 = child % childarray(1,1,2,nb,ndir)
  371. i2 = child % childarray(1,2,2,nb,ndir)
  372. j1 = child % childarray(2,1,2,nb,ndir)
  373. j2 = child % childarray(2,2,2,nb,ndir)
  374. k1 = child % childarray(3,1,2,nb,ndir)
  375. k2 = child % childarray(3,2,2,nb,ndir)
  376. l1 = child % childarray(4,1,2,nb,ndir)
  377. l2 = child % childarray(4,2,2,nb,ndir)
  378. m1 = child % childarray(5,1,2,nb,ndir)
  379. m2 = child % childarray(5,2,2,nb,ndir)
  380. n1 = child % childarray(6,1,2,nb,ndir)
  381. n2 = child % childarray(6,2,2,nb,ndir)
  382. call procname(parray6(i1:i2,j1:j2,k1:k2,l1:l2,m1:m2,n1:n2), &
  383. i1,i2,j1,j2,k1,k2,l1,l2,m1,m2,n1,n2, .FALSE.,coords(nb),ndir)
  384. end select
  385. endif
  386. enddo
  387. enddo
  388. !---------------------------------------------------------------------------------------------------
  389. end subroutine Agrif_Correctnd
  390. !===================================================================================================
  391. !
  392. !===================================================================================================
  393. ! subroutine saveAfterInterp
  394. !
  395. !> saves the values of the grid variable on the fine grid after the space interpolation
  396. !---------------------------------------------------------------------------------------------------
  397. subroutine saveAfterInterp ( child_var, bounds, kindex, newsize, nbdim )
  398. !---------------------------------------------------------------------------------------------------
  399. TYPE (Agrif_Variable), INTENT(inout) :: child_var !< The fine grid variable
  400. INTEGER, DIMENSION(nbdim,2), INTENT(in) :: bounds
  401. INTEGER, INTENT(inout) :: kindex !< Index indicating where this safeguard
  402. !< is done on the fine grid
  403. INTEGER, INTENT(in) :: newsize
  404. INTEGER, INTENT(in) :: nbdim
  405. !
  406. INTEGER :: ir,jr,kr,lr,mr,nr
  407. !
  408. ! Allocation of the array oldvalues2d
  409. !
  410. if (newsize .LE. 0) return
  411. !
  412. Call Agrif_Checksize(child_var,kindex+newsize)
  413. if (child_var % interpIndex /= Agrif_Curgrid % parent % ngridstep ) then
  414. child_var % oldvalues2d(1,kindex:kindex+newsize-1) = &
  415. child_var % oldvalues2d(2,kindex:kindex+newsize-1)
  416. endif
  417. SELECT CASE (nbdim)
  418. CASE (1)
  419. !CDIR ALTCODE
  420. do ir = bounds(1,1), bounds(1,2)
  421. child_var % oldvalues2d(2,kindex) = parray1(ir)
  422. kindex = kindex + 1
  423. enddo
  424. !
  425. CASE (2)
  426. do jr = bounds(2,1),bounds(2,2)
  427. !CDIR ALTCODE
  428. do ir = bounds(1,1),bounds(1,2)
  429. child_var % oldvalues2d(2,kindex) = parray2(ir,jr)
  430. kindex = kindex + 1
  431. enddo
  432. enddo
  433. !
  434. CASE (3)
  435. do kr = bounds(3,1),bounds(3,2)
  436. do jr = bounds(2,1),bounds(2,2)
  437. !CDIR ALTCODE
  438. do ir = bounds(1,1),bounds(1,2)
  439. child_var % oldvalues2d(2,kindex) = parray3(ir,jr,kr)
  440. kindex = kindex + 1
  441. enddo
  442. enddo
  443. enddo
  444. !
  445. CASE (4)
  446. do lr = bounds(4,1),bounds(4,2)
  447. do kr = bounds(3,1),bounds(3,2)
  448. do jr = bounds(2,1),bounds(2,2)
  449. !CDIR ALTCODE
  450. do ir = bounds(1,1),bounds(1,2)
  451. child_var % oldvalues2d(2,kindex) = parray4(ir,jr,kr,lr)
  452. kindex = kindex + 1
  453. enddo
  454. enddo
  455. enddo
  456. enddo
  457. !
  458. CASE (5)
  459. do mr = bounds(5,1),bounds(5,2)
  460. do lr = bounds(4,1),bounds(4,2)
  461. do kr = bounds(3,1),bounds(3,2)
  462. do jr = bounds(2,1),bounds(2,2)
  463. !CDIR ALTCODE
  464. do ir = bounds(1,1),bounds(1,2)
  465. child_var % oldvalues2d(2,kindex) = parray5(ir,jr,kr,lr,mr)
  466. kindex = kindex + 1
  467. enddo
  468. enddo
  469. enddo
  470. enddo
  471. enddo
  472. !
  473. CASE (6)
  474. do nr = bounds(6,1),bounds(6,2)
  475. do mr = bounds(5,1),bounds(5,2)
  476. do lr = bounds(4,1),bounds(4,2)
  477. do kr = bounds(3,1),bounds(3,2)
  478. do jr = bounds(2,1),bounds(2,2)
  479. !CDIR ALTCODE
  480. do ir = bounds(1,1),bounds(1,2)
  481. child_var % oldvalues2d(2,kindex) = parray6(ir,jr,kr,lr,mr,nr)
  482. kindex = kindex + 1
  483. enddo
  484. enddo
  485. enddo
  486. enddo
  487. enddo
  488. enddo
  489. END SELECT
  490. !---------------------------------------------------------------------------------------------------
  491. end subroutine saveAfterInterp
  492. !===================================================================================================
  493. !
  494. !===================================================================================================
  495. ! subroutine timeInterpolation
  496. !
  497. !> subroutine for a linear time interpolation on the child grid
  498. !---------------------------------------------------------------------------------------------------
  499. subroutine timeInterpolation ( child_var, bounds, kindex, c1t, c2t, nbdim )
  500. !---------------------------------------------------------------------------------------------------
  501. TYPE (Agrif_Variable) :: child_var !< The fine grid variable
  502. INTEGER, DIMENSION(nbdim,2) :: bounds
  503. INTEGER :: kindex !< Index indicating the values of the fine grid got
  504. !< before and after the space interpolation and
  505. !< used for the time interpolation
  506. REAL :: c1t, c2t !< Coefficients for the time interpolation (c2t=1-c1t)
  507. INTEGER :: nbdim
  508. !
  509. INTEGER :: ir,jr,kr,lr,mr,nr
  510. !
  511. SELECT CASE (nbdim)
  512. CASE (1)
  513. !CDIR ALTCODE
  514. do ir = bounds(1,1),bounds(1,2)
  515. parray1(ir) = c2t*child_var % oldvalues2d(1,kindex) + &
  516. c1t*child_var % oldvalues2d(2,kindex)
  517. kindex = kindex + 1
  518. enddo
  519. !
  520. CASE (2)
  521. do jr = bounds(2,1),bounds(2,2)
  522. !CDIR ALTCODE
  523. do ir = bounds(1,1),bounds(1,2)
  524. parray2(ir,jr) = c2t*child_var % oldvalues2d(1,kindex) + &
  525. c1t*child_var % oldvalues2d(2,kindex)
  526. kindex = kindex + 1
  527. enddo
  528. enddo
  529. !
  530. CASE (3)
  531. do kr = bounds(3,1),bounds(3,2)
  532. do jr = bounds(2,1),bounds(2,2)
  533. !CDIR ALTCODE
  534. do ir = bounds(1,1),bounds(1,2)
  535. parray3(ir,jr,kr) = c2t*child_var % oldvalues2d(1,kindex) + &
  536. c1t*child_var % oldvalues2d(2,kindex)
  537. kindex = kindex + 1
  538. enddo
  539. enddo
  540. enddo
  541. !
  542. CASE (4)
  543. do lr = bounds(4,1),bounds(4,2)
  544. do kr = bounds(3,1),bounds(3,2)
  545. do jr = bounds(2,1),bounds(2,2)
  546. !CDIR ALTCODE
  547. do ir = bounds(1,1),bounds(1,2)
  548. parray4(ir,jr,kr,lr) = c2t*child_var % oldvalues2d(1,kindex) + &
  549. c1t*child_var % oldvalues2d(2,kindex)
  550. kindex = kindex + 1
  551. enddo
  552. enddo
  553. enddo
  554. enddo
  555. !
  556. CASE (5)
  557. do mr=bounds(5,1),bounds(5,2)
  558. do lr=bounds(4,1),bounds(4,2)
  559. do kr=bounds(3,1),bounds(3,2)
  560. do jr=bounds(2,1),bounds(2,2)
  561. !CDIR ALTCODE
  562. do ir=bounds(1,1),bounds(1,2)
  563. parray5(ir,jr,kr,lr,mr) = c2t*child_var % oldvalues2d(1,kindex) + &
  564. c1t*child_var % oldvalues2d(2,kindex)
  565. kindex = kindex + 1
  566. enddo
  567. enddo
  568. enddo
  569. enddo
  570. enddo
  571. !
  572. CASE (6)
  573. do nr=bounds(6,1),bounds(6,2)
  574. do mr=bounds(5,1),bounds(5,2)
  575. do lr=bounds(4,1),bounds(4,2)
  576. do kr=bounds(3,1),bounds(3,2)
  577. do jr=bounds(2,1),bounds(2,2)
  578. !CDIR ALTCODE
  579. do ir=bounds(1,1),bounds(1,2)
  580. parray6(ir,jr,kr,lr,mr,nr) = c2t*child_var % oldvalues2d(1,kindex) + &
  581. c1t*child_var % oldvalues2d(2,kindex)
  582. kindex = kindex + 1
  583. enddo
  584. enddo
  585. enddo
  586. enddo
  587. enddo
  588. enddo
  589. END SELECT
  590. !---------------------------------------------------------------------------------------------------
  591. end subroutine timeInterpolation
  592. !===================================================================================================
  593. !
  594. !===================================================================================================
  595. ! subroutine Agrif_Checksize
  596. !
  597. !> subroutine used in the saveAfterInterp procedure to allocate the oldvalues2d array
  598. !---------------------------------------------------------------------------------------------------
  599. subroutine Agrif_Checksize ( child_var, newsize )
  600. !---------------------------------------------------------------------------------------------------
  601. TYPE (Agrif_Variable), INTENT(inout) :: child_var !< The fine grid variable
  602. INTEGER , INTENT(in) :: newsize !< Size of the domains where the boundary
  603. !< conditions are calculated
  604. !
  605. REAL, DIMENSION(:,:), Allocatable :: tempoldvalues ! Temporary array
  606. !
  607. if (.NOT. associated(child_var % oldvalues2d)) then
  608. !
  609. allocate(child_var % oldvalues2d(2,newsize))
  610. child_var % oldvalues2d = 0.
  611. !
  612. else
  613. !
  614. if (SIZE(child_var % oldvalues2d,2) < newsize) then
  615. !
  616. allocate(tempoldvalues(2,SIZE(child_var % oldvalues2d,2)))
  617. tempoldvalues = child_var % oldvalues2d
  618. deallocate(child_var % oldvalues2d)
  619. allocate( child_var % oldvalues2d(2,newsize))
  620. child_var % oldvalues2d = 0.
  621. child_var % oldvalues2d(:,1:SIZE(tempoldvalues,2)) = tempoldvalues(:,:)
  622. deallocate(tempoldvalues)
  623. !
  624. endif
  625. !
  626. endif
  627. !---------------------------------------------------------------------------------------------------
  628. end subroutine Agrif_Checksize
  629. !===================================================================================================
  630. !
  631. end module Agrif_Boundary