modcurgridfunctions.F90 36 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765
  1. !
  2. ! $Id: modcurgridfunctions.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 to define some procedures concerning the current grid
  24. !
  25. module Agrif_CurgridFunctions
  26. !
  27. use Agrif_Init
  28. !
  29. implicit none
  30. !
  31. contains
  32. !
  33. !===================================================================================================
  34. ! function Agrif_rel_dt
  35. !
  36. !> Returns the time step of the current grid, relatively to the root grid (for which dt=1.).
  37. !---------------------------------------------------------------------------------------------------
  38. function Agrif_rel_dt ( ) result( rel_dt )
  39. !---------------------------------------------------------------------------------------------------
  40. integer :: i
  41. real :: rel_dt
  42. !
  43. rel_dt = 1.
  44. !
  45. do i = 1,Agrif_Probdim
  46. rel_dt = min(rel_dt, Agrif_Curgrid % Agrif_dt(i))
  47. enddo
  48. !---------------------------------------------------------------------------------------------------
  49. end function Agrif_rel_dt
  50. !===================================================================================================
  51. !
  52. !===================================================================================================
  53. ! function Agrif_rel_idt
  54. !
  55. !> Returns the time refinement factor of the current grid, relatively to the root grid (for which idt=1).
  56. !---------------------------------------------------------------------------------------------------
  57. function Agrif_rel_idt ( ) result( rel_idt )
  58. !---------------------------------------------------------------------------------------------------
  59. integer :: rel_idt
  60. !
  61. rel_idt = nint(1./Agrif_rel_dt())
  62. !---------------------------------------------------------------------------------------------------
  63. end function Agrif_rel_idt
  64. !===================================================================================================
  65. !
  66. !===================================================================================================
  67. ! function Agrif_IRhot
  68. !
  69. !> Returns the time refinement factor of the current grid.
  70. !---------------------------------------------------------------------------------------------------
  71. function Agrif_IRhot ( ) result( irhot )
  72. !---------------------------------------------------------------------------------------------------
  73. integer :: i, irhot
  74. !
  75. irhot = 1
  76. !
  77. do i = 1,Agrif_Probdim
  78. irhot = max(irhot, Agrif_Curgrid % timeref(i))
  79. enddo
  80. !---------------------------------------------------------------------------------------------------
  81. end function Agrif_IRhot
  82. !===================================================================================================
  83. !
  84. !===================================================================================================
  85. ! function Agrif_Rhot
  86. !
  87. !> Returns the time refinement factor of the current grid.
  88. !---------------------------------------------------------------------------------------------------
  89. function Agrif_Rhot ( ) result( rhot )
  90. !---------------------------------------------------------------------------------------------------
  91. real :: rhot
  92. !
  93. rhot = float(Agrif_IRhot())
  94. !---------------------------------------------------------------------------------------------------
  95. end function Agrif_Rhot
  96. !===================================================================================================
  97. !
  98. !===================================================================================================
  99. ! function Agrif_Parent_IRhot
  100. !
  101. !> Returns the time refinement factor of the parent of the current grid.
  102. !---------------------------------------------------------------------------------------------------
  103. function Agrif_Parent_IRhot ( ) result( irhot )
  104. !---------------------------------------------------------------------------------------------------
  105. integer :: i, irhot
  106. !
  107. irhot = 1
  108. !
  109. do i = 1,Agrif_Probdim
  110. irhot = max(irhot, Agrif_Curgrid % parent % timeref(i))
  111. enddo
  112. !---------------------------------------------------------------------------------------------------
  113. end function Agrif_Parent_IRhot
  114. !===================================================================================================
  115. !
  116. !===================================================================================================
  117. ! function Agrif_Parent_Rhot
  118. !
  119. !> Returns the time refinement factor of the parent of the current grid.
  120. !---------------------------------------------------------------------------------------------------
  121. function Agrif_Parent_Rhot ( ) result( rhot )
  122. !---------------------------------------------------------------------------------------------------
  123. real :: rhot
  124. !
  125. rhot = float(Agrif_Parent_IRhot())
  126. !---------------------------------------------------------------------------------------------------
  127. end function Agrif_Parent_Rhot
  128. !===================================================================================================
  129. !
  130. !===================================================================================================
  131. ! function Agrif_Nbstepint
  132. !
  133. !> function for the calculation of the coefficients used for the time interpolation
  134. !! (module #Agrif_Boundary).
  135. !---------------------------------------------------------------------------------------------------
  136. function Agrif_Nbstepint ( )
  137. !---------------------------------------------------------------------------------------------------
  138. integer :: Agrif_nbstepint ! result
  139. !
  140. Agrif_nbstepint = mod(Agrif_Curgrid % ngridstep, Agrif_iRhot())
  141. !---------------------------------------------------------------------------------------------------
  142. end function Agrif_Nbstepint
  143. !===================================================================================================
  144. !
  145. !===================================================================================================
  146. ! function Agrif_Parent_Nbstepint
  147. !
  148. !> function for the calculation of the coefficients used for the time interpolation
  149. !! (module #Agrif_Boundary).
  150. !---------------------------------------------------------------------------------------------------
  151. function Agrif_Parent_Nbstepint ( )
  152. !---------------------------------------------------------------------------------------------------
  153. integer :: Agrif_Parent_Nbstepint ! result
  154. !
  155. Agrif_Parent_Nbstepint = mod(Agrif_Curgrid % parent % ngridstep, int(Agrif_Parent_Rhot()))
  156. !---------------------------------------------------------------------------------------------------
  157. end function Agrif_Parent_Nbstepint
  158. !===================================================================================================
  159. !
  160. !===================================================================================================
  161. ! subroutine Agrif_InterpNearBorderX
  162. !
  163. !> Allows to interpolate (in the x direction) on a near border of the current grid if this one
  164. !! has a common border with the root coarse grid.
  165. !---------------------------------------------------------------------------------------------------
  166. subroutine Agrif_InterpNearBorderX ( )
  167. !---------------------------------------------------------------------------------------------------
  168. Agrif_Curgrid % NearRootBorder(1) = .FALSE.
  169. !---------------------------------------------------------------------------------------------------
  170. end subroutine Agrif_InterpNearBorderX
  171. !===================================================================================================
  172. !
  173. !===================================================================================================
  174. ! subroutine Agrif_InterpDistantBorderX
  175. !
  176. !> Allows to interpolate (in the x direction) on a distant border of the current grid if this one
  177. !! has a common border with the root coarse grid.
  178. !---------------------------------------------------------------------------------------------------
  179. subroutine Agrif_InterpDistantBorderX ( )
  180. !---------------------------------------------------------------------------------------------------
  181. Agrif_Curgrid % DistantRootBorder(1) = .FALSE.
  182. !---------------------------------------------------------------------------------------------------
  183. end subroutine Agrif_InterpDistantBorderX
  184. !===================================================================================================
  185. !
  186. !===================================================================================================
  187. ! subroutine Agrif_InterpNearBorderY
  188. !
  189. !> Allows to interpolate (in the y direction) on a near border of the current grid if this one
  190. !! has a common border with the root coarse grid.
  191. !---------------------------------------------------------------------------------------------------
  192. subroutine Agrif_InterpNearBorderY ( )
  193. !---------------------------------------------------------------------------------------------------
  194. Agrif_Curgrid % NearRootBorder(2) = .FALSE.
  195. !---------------------------------------------------------------------------------------------------
  196. end subroutine Agrif_InterpNearBorderY
  197. !===================================================================================================
  198. !
  199. !===================================================================================================
  200. ! subroutine Agrif_InterpDistantBorderY
  201. !
  202. !> Allows to interpolate (in the y direction) on a distant border of the current grid if this one
  203. !! has a common border with the root coarse grid.
  204. !---------------------------------------------------------------------------------------------------
  205. subroutine Agrif_InterpDistantBorderY ( )
  206. !---------------------------------------------------------------------------------------------------
  207. Agrif_Curgrid % DistantRootBorder(2) = .FALSE.
  208. !---------------------------------------------------------------------------------------------------
  209. end subroutine Agrif_InterpDistantBorderY
  210. !===================================================================================================
  211. !
  212. !===================================================================================================
  213. ! subroutine Agrif_InterpNearBorderZ
  214. !
  215. !> Allows to interpolate (in the z direction) on a near border of the current grid if this one
  216. !! has a common border with the root coarse grid.
  217. !---------------------------------------------------------------------------------------------------
  218. subroutine Agrif_InterpNearBorderZ ( )
  219. !---------------------------------------------------------------------------------------------------
  220. Agrif_Curgrid % NearRootBorder(3) = .FALSE.
  221. !---------------------------------------------------------------------------------------------------
  222. end subroutine Agrif_InterpNearBorderZ
  223. !===================================================================================================
  224. !
  225. !===================================================================================================
  226. ! subroutine Agrif_InterpDistantBorderZ
  227. !
  228. !> Allows to interpolate (in the z direction) on a distant border of the current grid if this one
  229. !! has a common border with the root coarse grid.
  230. !---------------------------------------------------------------------------------------------------
  231. subroutine Agrif_InterpDistantBorderZ()
  232. !---------------------------------------------------------------------------------------------------
  233. Agrif_Curgrid % DistantRootBorder(3) = .FALSE.
  234. !---------------------------------------------------------------------------------------------------
  235. end subroutine Agrif_InterpDistantBorderZ
  236. !===================================================================================================
  237. !
  238. !===================================================================================================
  239. ! function Agrif_Parent_Nb_Step
  240. !
  241. !> Returns the number of time steps of the parent of the current grid.
  242. !---------------------------------------------------------------------------------------------------
  243. function Agrif_Parent_Nb_Step ( )
  244. !---------------------------------------------------------------------------------------------------
  245. integer :: Agrif_Parent_Nb_Step ! Result
  246. !
  247. if (Agrif_Root()) then
  248. Agrif_Parent_Nb_Step = -1
  249. else
  250. Agrif_Parent_Nb_Step = Agrif_Curgrid % parent % ngridstep
  251. endif
  252. !---------------------------------------------------------------------------------------------------
  253. end function Agrif_Parent_Nb_Step
  254. !===================================================================================================
  255. !
  256. !===================================================================================================
  257. ! function Agrif_Root
  258. !
  259. !> Indicates if the current grid is or not the root grid.
  260. !---------------------------------------------------------------------------------------------------
  261. function Agrif_Root ( )
  262. !---------------------------------------------------------------------------------------------------
  263. logical :: Agrif_Root ! Result
  264. !
  265. Agrif_Root = (Agrif_Curgrid % fixedrank == 0)
  266. !---------------------------------------------------------------------------------------------------
  267. end function Agrif_Root
  268. !===================================================================================================
  269. !
  270. !===================================================================================================
  271. ! function Agrif_GrandMother
  272. !
  273. !> Indicates if the current grid is or not the root grid.
  274. !---------------------------------------------------------------------------------------------------
  275. function Agrif_GrandMother ( )
  276. !---------------------------------------------------------------------------------------------------
  277. logical :: Agrif_GrandMother ! Result
  278. !
  279. Agrif_GrandMother = Agrif_Curgrid % grand_mother_grid
  280. !---------------------------------------------------------------------------------------------------
  281. end function Agrif_GrandMother
  282. !===================================================================================================
  283. !
  284. !===================================================================================================
  285. ! function Agrif_Parent_Root
  286. !
  287. !> Indicates if the parent of the current grid is or not the root grid.
  288. !---------------------------------------------------------------------------------------------------
  289. function Agrif_Parent_Root ( )
  290. !---------------------------------------------------------------------------------------------------
  291. logical :: Agrif_Parent_Root ! Result
  292. !
  293. Agrif_Parent_Root = (Agrif_Curgrid % parent % fixedrank == 0)
  294. !---------------------------------------------------------------------------------------------------
  295. end function Agrif_Parent_Root
  296. !===================================================================================================
  297. !
  298. !===================================================================================================
  299. ! function Agrif_Fixed
  300. !
  301. !> Returns the number of the current grid.
  302. !---------------------------------------------------------------------------------------------------
  303. function Agrif_Fixed ( )
  304. !---------------------------------------------------------------------------------------------------
  305. integer :: Agrif_Fixed ! Result
  306. !
  307. if (Agrif_Curgrid % fixed) then
  308. Agrif_Fixed = Agrif_Curgrid % fixedrank
  309. else
  310. Agrif_Fixed = -1
  311. endif
  312. !---------------------------------------------------------------------------------------------------
  313. end function Agrif_Fixed
  314. !===================================================================================================
  315. !
  316. !===================================================================================================
  317. ! function Agrif_Parent_Fixed
  318. !
  319. !> Returns the number of the parent of the current grid.
  320. !---------------------------------------------------------------------------------------------------
  321. function Agrif_Parent_Fixed ( )
  322. !---------------------------------------------------------------------------------------------------
  323. integer :: Agrif_Parent_Fixed ! Result
  324. !
  325. if (Agrif_Curgrid % parent % fixed) then
  326. Agrif_Parent_Fixed = Agrif_Curgrid % parent % fixedrank
  327. else
  328. Agrif_Parent_Fixed = 0
  329. endif
  330. !---------------------------------------------------------------------------------------------------
  331. end function Agrif_Parent_Fixed
  332. !===================================================================================================
  333. !
  334. !===================================================================================================
  335. ! function Agrif_Is_Fixed
  336. !
  337. !> Returns .TRUE. if the current grid is fixed.
  338. !---------------------------------------------------------------------------------------------------
  339. function Agrif_Is_Fixed ( )
  340. !---------------------------------------------------------------------------------------------------
  341. logical :: Agrif_Is_Fixed ! Result
  342. !
  343. Agrif_Is_Fixed = Agrif_Curgrid % fixed
  344. !---------------------------------------------------------------------------------------------------
  345. end function Agrif_Is_Fixed
  346. !===================================================================================================
  347. !
  348. !===================================================================================================
  349. ! function Agrif_Parent_Is_Fixed
  350. !
  351. !> Returns .TRUE. if the parent of the current grid is fixed.
  352. !---------------------------------------------------------------------------------------------------
  353. function Agrif_Parent_Is_Fixed ( )
  354. !---------------------------------------------------------------------------------------------------
  355. logical :: Agrif_Parent_Is_Fixed ! Result
  356. !
  357. Agrif_Parent_Is_Fixed = Agrif_Curgrid % parent % fixed
  358. !---------------------------------------------------------------------------------------------------
  359. end function Agrif_Parent_Is_Fixed
  360. !===================================================================================================
  361. !
  362. !===================================================================================================
  363. ! function Agrif_CFixed
  364. !
  365. !> Returns the number of the current grid.
  366. !---------------------------------------------------------------------------------------------------
  367. function Agrif_CFixed ( )
  368. !---------------------------------------------------------------------------------------------------
  369. character(3) :: Agrif_CFixed ! Result
  370. !
  371. character(3) :: cfixed
  372. integer :: fixed
  373. !
  374. fixed = Agrif_Fixed()
  375. !
  376. if (fixed /= -1) then
  377. !
  378. if (fixed <= 9) then
  379. write(cfixed,'(i1)') fixed
  380. else
  381. write(cfixed,'(i2)') fixed
  382. endif
  383. !
  384. Agrif_CFixed = cfixed
  385. if (associated(agrif_curgrid,agrif_coarsegrid)) then
  386. Agrif_CFixed = 'gm'
  387. endif
  388. !
  389. else
  390. print*,'Call to Agrif_CFixed() on a moving grid'
  391. stop
  392. endif
  393. !---------------------------------------------------------------------------------------------------
  394. end function Agrif_CFixed
  395. !===================================================================================================
  396. !
  397. !===================================================================================================
  398. ! function Agrid_Parent_CFixed
  399. !
  400. !> Returns the number of the parent of the current grid.
  401. !---------------------------------------------------------------------------------------------------
  402. function Agrid_Parent_CFixed ( )
  403. !---------------------------------------------------------------------------------------------------
  404. character(3) :: Agrid_Parent_CFixed ! Result
  405. !
  406. character(3) :: cfixed
  407. integer :: fixed
  408. !
  409. fixed = Agrif_Parent_Fixed()
  410. !
  411. if(fixed /= -1) then
  412. !
  413. if (fixed <= 9) then
  414. write(cfixed,'(i1)')fixed
  415. else
  416. write(cfixed,'(i2)')fixed
  417. endif
  418. !
  419. Agrid_Parent_CFixed=cfixed
  420. !
  421. else
  422. print*,'Illegal call to Agrid_Parent_CFixed()'
  423. stop
  424. endif
  425. !---------------------------------------------------------------------------------------------------
  426. end function Agrid_Parent_CFixed
  427. !===================================================================================================
  428. !
  429. !===================================================================================================
  430. ! subroutine Agrif_ChildGrid_to_ParentGrid
  431. !
  432. !> Make the pointer #Agrif_Curgrid point on the parent grid of the current grid.
  433. !---------------------------------------------------------------------------------------------------
  434. subroutine Agrif_ChildGrid_to_ParentGrid ( )
  435. !---------------------------------------------------------------------------------------------------
  436. Agrif_Curgrid % parent % save_grid => Agrif_Curgrid
  437. call Agrif_Instance(Agrif_Curgrid%parent)
  438. !---------------------------------------------------------------------------------------------------
  439. end subroutine Agrif_ChildGrid_to_ParentGrid
  440. !===================================================================================================
  441. !
  442. !===================================================================================================
  443. ! subroutine Agrif_ParentGrid_to_ChildGrid
  444. !
  445. !> Make the pointer #Agrif_Curgrid point on the child grid after having called the
  446. !! #Agrif_ChildGrid_to_ParentGrid subroutine.
  447. !---------------------------------------------------------------------------------------------------
  448. subroutine Agrif_ParentGrid_to_ChildGrid ( )
  449. !---------------------------------------------------------------------------------------------------
  450. call Agrif_Instance(Agrif_Curgrid%save_grid)
  451. !---------------------------------------------------------------------------------------------------
  452. end subroutine Agrif_ParentGrid_to_ChildGrid
  453. !===================================================================================================
  454. !
  455. !===================================================================================================
  456. ! function Agrif_Get_Unit
  457. !
  458. !> Returns a unit not connected to any file.
  459. !---------------------------------------------------------------------------------------------------
  460. function Agrif_Get_Unit ( )
  461. !---------------------------------------------------------------------------------------------------
  462. integer :: Agrif_Get_Unit ! Result
  463. !
  464. integer :: n
  465. logical :: op
  466. !
  467. integer :: nunit
  468. integer :: iii, out, iiimax
  469. logical :: bexist
  470. integer,dimension(1:1000) :: forbiddenunit
  471. !
  472. ! Load forbidden Unit if the file Agrif_forbidenUnit exist
  473. !
  474. INQUIRE(file='Agrif_forbiddenUnit.txt', exist=bexist)
  475. !
  476. if (.not. bexist) then
  477. ! File Agrif_forbiddenUnit.txt not found
  478. else
  479. nunit = 777
  480. OPEN(nunit,file='Agrif_forbiddenUnit.txt', form='formatted', status="old")
  481. iii = 1
  482. do while ( .TRUE. )
  483. READ(nunit,*, end=99) forbiddenunit(iii)
  484. iii = iii + 1
  485. enddo
  486. 99 continue
  487. iiimax = iii
  488. close(nunit)
  489. endif
  490. !
  491. do n = 7,1000
  492. !
  493. INQUIRE(Unit=n,Opened=op)
  494. !
  495. out = 0
  496. if ( bexist .AND. (.NOT.op) ) then
  497. do iii = 1,iiimax
  498. if ( n == forbiddenunit(iii) ) out = 1
  499. enddo
  500. endif
  501. !
  502. if ( (.NOT.op) .AND. (out == 0) ) exit
  503. !
  504. enddo
  505. !
  506. Agrif_Get_Unit = n
  507. !---------------------------------------------------------------------------------------------------
  508. end function Agrif_Get_Unit
  509. !===================================================================================================
  510. !
  511. !===================================================================================================
  512. ! subroutine Agrif_Set_Extra_Boundary_Cells
  513. !---------------------------------------------------------------------------------------------------
  514. subroutine Agrif_Set_Extra_Boundary_Cells ( nb_extra_cells )
  515. !---------------------------------------------------------------------------------------------------
  516. integer, intent(in) :: nb_extra_cells
  517. !
  518. Agrif_Extra_Boundary_Cells = nb_extra_cells
  519. !---------------------------------------------------------------------------------------------------
  520. end subroutine Agrif_Set_Extra_Boundary_Cells
  521. !===================================================================================================
  522. !
  523. !===================================================================================================
  524. ! subroutine Agrif_Set_Efficiency
  525. !---------------------------------------------------------------------------------------------------
  526. subroutine Agrif_Set_Efficiency ( eff )
  527. !---------------------------------------------------------------------------------------------------
  528. real, intent(in) :: eff
  529. !
  530. if ( (eff < 0.) .OR. (eff > 1) ) then
  531. write(*,*) 'Error Efficiency should be between 0 and 1'
  532. stop
  533. else
  534. Agrif_Efficiency = eff
  535. endif
  536. !---------------------------------------------------------------------------------------------------
  537. end subroutine Agrif_Set_Efficiency
  538. !===================================================================================================
  539. !
  540. !===================================================================================================
  541. ! subroutine Agrif_Set_Regridding
  542. !---------------------------------------------------------------------------------------------------
  543. subroutine Agrif_Set_Regridding ( regfreq )
  544. !---------------------------------------------------------------------------------------------------
  545. integer, intent(in) :: regfreq
  546. !
  547. if (regfreq < 0) then
  548. write(*,*) 'Regridding frequency should be positive'
  549. stop
  550. else
  551. Agrif_Regridding = regfreq
  552. endif
  553. !---------------------------------------------------------------------------------------------------
  554. end subroutine Agrif_Set_Regridding
  555. !===================================================================================================
  556. !
  557. !===================================================================================================
  558. ! subroutine Agrif_Set_coeffref_x
  559. !---------------------------------------------------------------------------------------------------
  560. subroutine Agrif_Set_coeffref_x ( coeffref )
  561. !---------------------------------------------------------------------------------------------------
  562. integer, intent(in) :: coeffref
  563. if (coeffref < 0) then
  564. write(*,*) 'Coefficient of raffinement should be positive'
  565. stop
  566. else
  567. Agrif_coeffref(1) = coeffref
  568. endif
  569. !---------------------------------------------------------------------------------------------------
  570. end subroutine Agrif_Set_coeffref_x
  571. !===================================================================================================
  572. !
  573. !===================================================================================================
  574. ! subroutine Agrif_Set_coeffref_y
  575. !---------------------------------------------------------------------------------------------------
  576. subroutine Agrif_Set_coeffref_y ( coeffref )
  577. !---------------------------------------------------------------------------------------------------
  578. integer, intent(in) :: coeffref
  579. if (coeffref < 0) then
  580. write(*,*) 'Coefficient of raffinement should be positive'
  581. stop
  582. else
  583. Agrif_coeffref(2) = coeffref
  584. endif
  585. !---------------------------------------------------------------------------------------------------
  586. end subroutine Agrif_Set_coeffref_y
  587. !===================================================================================================
  588. !
  589. !===================================================================================================
  590. ! subroutine Agrif_Set_coeffref_z
  591. !---------------------------------------------------------------------------------------------------
  592. subroutine Agrif_Set_coeffref_z ( coeffref )
  593. !---------------------------------------------------------------------------------------------------
  594. integer, intent(in) :: coeffref
  595. !
  596. if (coeffref < 0) then
  597. write(*,*) 'Coefficient of raffinement should be positive'
  598. stop
  599. else
  600. Agrif_coeffref(3) = coeffref
  601. endif
  602. !---------------------------------------------------------------------------------------------------
  603. end subroutine Agrif_Set_coeffref_z
  604. !===================================================================================================
  605. !
  606. !===================================================================================================
  607. ! subroutine Agrif_Set_coeffreft_x
  608. !---------------------------------------------------------------------------------------------------
  609. subroutine Agrif_Set_coeffreft_x ( coeffref )
  610. !---------------------------------------------------------------------------------------------------
  611. integer, intent(in) :: coeffref
  612. if (coeffref < 0) then
  613. write(*,*) 'Coefficient of time raffinement should be positive'
  614. stop
  615. else
  616. Agrif_coeffreft(1) = coeffref
  617. endif
  618. !---------------------------------------------------------------------------------------------------
  619. end subroutine Agrif_Set_coeffreft_x
  620. !===================================================================================================
  621. !
  622. !===================================================================================================
  623. ! subroutine Agrif_Set_coeffreft_y
  624. !---------------------------------------------------------------------------------------------------
  625. subroutine Agrif_Set_coeffreft_y ( coeffref )
  626. !---------------------------------------------------------------------------------------------------
  627. integer, intent(in) :: coeffref
  628. !
  629. if (coeffref < 0) then
  630. write(*,*) 'Coefficient of time raffinement should be positive'
  631. stop
  632. else
  633. Agrif_coeffreft(2) = coeffref
  634. endif
  635. !---------------------------------------------------------------------------------------------------
  636. end subroutine Agrif_Set_coeffreft_y
  637. !===================================================================================================
  638. !
  639. !===================================================================================================
  640. ! subroutine Agrif_Set_coeffreft_z
  641. !---------------------------------------------------------------------------------------------------
  642. subroutine Agrif_Set_coeffreft_z ( coeffref )
  643. !---------------------------------------------------------------------------------------------------
  644. integer, intent(in) :: coeffref
  645. if (coeffref < 0) then
  646. write(*,*)'Coefficient of time raffinement should be positive'
  647. stop
  648. else
  649. Agrif_coeffreft(3) = coeffref
  650. endif
  651. !---------------------------------------------------------------------------------------------------
  652. end subroutine Agrif_Set_coeffreft_z
  653. !===================================================================================================
  654. !
  655. !===================================================================================================
  656. ! subroutine Agrif_Set_Minwidth
  657. !---------------------------------------------------------------------------------------------------
  658. subroutine Agrif_Set_Minwidth ( coefminwidth )
  659. !---------------------------------------------------------------------------------------------------
  660. integer, intent(in) :: coefminwidth
  661. !
  662. if (coefminwidth < 0) then
  663. write(*,*)'Coefficient of Minwidth should be positive'
  664. stop
  665. else
  666. Agrif_Minwidth = coefminwidth
  667. endif
  668. !---------------------------------------------------------------------------------------------------
  669. end subroutine Agrif_Set_Minwidth
  670. !===================================================================================================
  671. !
  672. !===================================================================================================
  673. ! subroutine Agrif_Set_Rafmax
  674. !---------------------------------------------------------------------------------------------------
  675. subroutine Agrif_Set_Rafmax ( coefrafmax )
  676. !---------------------------------------------------------------------------------------------------
  677. integer, intent(in) :: coefrafmax
  678. !
  679. integer :: i
  680. real :: res
  681. !
  682. if (coefrafmax < 0) then
  683. write(*,*)'Coefficient of should be positive'
  684. stop
  685. else
  686. res = 1.
  687. do i = 1,coefrafmax-1
  688. res = res * FLOAT(Agrif_coeffref(1))
  689. enddo
  690. if ( res == 0 ) res = 1
  691. Agrif_Mind(1) = 1. / res
  692. !
  693. res = 1.
  694. do i = 1,coefrafmax-1
  695. res = res * FLOAT(Agrif_coeffref(2))
  696. enddo
  697. if ( res == 0 ) res = 1
  698. Agrif_Mind(2) = 1. / res
  699. !
  700. res = 1.
  701. do i = 1,coefrafmax-1
  702. res = res * FLOAT(Agrif_coeffref(3))
  703. enddo
  704. if ( res == 0 ) res = 1
  705. Agrif_Mind(3) = 1. / res
  706. !
  707. endif
  708. !---------------------------------------------------------------------------------------------------
  709. end subroutine Agrif_Set_Rafmax
  710. !===================================================================================================
  711. !
  712. !===================================================================================================
  713. ! subroutine Agrif_Set_MaskMaxSearch
  714. !---------------------------------------------------------------------------------------------------
  715. subroutine Agrif_Set_MaskMaxSearch ( mymaxsearch )
  716. !---------------------------------------------------------------------------------------------------
  717. integer, intent(in) :: mymaxsearch
  718. !
  719. MaxSearch = mymaxsearch
  720. !---------------------------------------------------------------------------------------------------
  721. end subroutine Agrif_Set_MaskMaxSearch
  722. !===================================================================================================
  723. !
  724. !===================================================================================================
  725. ! function Agrif_Level
  726. !---------------------------------------------------------------------------------------------------
  727. function Agrif_Level ( )
  728. !---------------------------------------------------------------------------------------------------
  729. integer :: Agrif_Level ! Result
  730. !
  731. Agrif_Level = Agrif_Curgrid % level
  732. !---------------------------------------------------------------------------------------------------
  733. end function Agrif_Level
  734. !===================================================================================================
  735. !
  736. !===================================================================================================
  737. ! function Agrif_MaxLevel
  738. !---------------------------------------------------------------------------------------------------
  739. function Agrif_MaxLevel ( )
  740. !---------------------------------------------------------------------------------------------------
  741. integer :: Agrif_MaxLevel ! Result
  742. !
  743. Agrif_MaxLevel = Agrif_MaxLevelLoc
  744. !---------------------------------------------------------------------------------------------------
  745. end function Agrif_MaxLevel
  746. !===================================================================================================
  747. !
  748. !===================================================================================================
  749. ! function Agrif_GridAllocation_is_done
  750. !---------------------------------------------------------------------------------------------------
  751. function Agrif_GridAllocation_is_done ( ) result(isdone)
  752. !---------------------------------------------------------------------------------------------------
  753. logical :: isdone
  754. !
  755. isdone = Agrif_Curgrid % allocation_is_done
  756. !---------------------------------------------------------------------------------------------------
  757. end function Agrif_GridAllocation_is_done
  758. !===================================================================================================
  759. !
  760. end module Agrif_CurgridFunctions