1 | ! |
---|
2 | ! $Id$ |
---|
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 | |
---|
32 | interface Agrif_Parent |
---|
33 | module procedure Agrif_Parent_Real_4, & |
---|
34 | Agrif_Parent_Real_8, & |
---|
35 | Agrif_Parent_Integer, & |
---|
36 | Agrif_Parent_Character, & |
---|
37 | Agrif_Parent_Logical |
---|
38 | end interface |
---|
39 | |
---|
40 | contains |
---|
41 | ! |
---|
42 | !=================================================================================================== |
---|
43 | ! function Agrif_rel_dt |
---|
44 | ! |
---|
45 | !> Returns the time step of the current grid, relatively to the root grid (for which dt=1.). |
---|
46 | !--------------------------------------------------------------------------------------------------- |
---|
47 | function Agrif_rel_dt ( ) result( rel_dt ) |
---|
48 | !--------------------------------------------------------------------------------------------------- |
---|
49 | integer :: i |
---|
50 | real :: rel_dt |
---|
51 | ! |
---|
52 | rel_dt = 1. |
---|
53 | ! |
---|
54 | do i = 1,Agrif_Probdim |
---|
55 | rel_dt = min(rel_dt, Agrif_Curgrid % Agrif_dt(i)) |
---|
56 | enddo |
---|
57 | !--------------------------------------------------------------------------------------------------- |
---|
58 | end function Agrif_rel_dt |
---|
59 | !=================================================================================================== |
---|
60 | ! |
---|
61 | !=================================================================================================== |
---|
62 | ! function Agrif_rel_idt |
---|
63 | ! |
---|
64 | !> Returns the time refinement factor of the current grid, relatively to the root grid (for which idt=1). |
---|
65 | !--------------------------------------------------------------------------------------------------- |
---|
66 | function Agrif_rel_idt ( ) result( rel_idt ) |
---|
67 | !--------------------------------------------------------------------------------------------------- |
---|
68 | integer :: rel_idt |
---|
69 | ! |
---|
70 | rel_idt = nint(1./Agrif_rel_dt()) |
---|
71 | !--------------------------------------------------------------------------------------------------- |
---|
72 | end function Agrif_rel_idt |
---|
73 | !=================================================================================================== |
---|
74 | ! |
---|
75 | !=================================================================================================== |
---|
76 | ! function Agrif_IRhot |
---|
77 | ! |
---|
78 | !> Returns the time refinement factor of the current grid. |
---|
79 | !--------------------------------------------------------------------------------------------------- |
---|
80 | function Agrif_IRhot ( ) result( irhot ) |
---|
81 | !--------------------------------------------------------------------------------------------------- |
---|
82 | integer :: i, irhot |
---|
83 | ! |
---|
84 | irhot = 1 |
---|
85 | ! |
---|
86 | do i = 1,Agrif_Probdim |
---|
87 | irhot = max(irhot, Agrif_Curgrid % timeref(i)) |
---|
88 | enddo |
---|
89 | !--------------------------------------------------------------------------------------------------- |
---|
90 | end function Agrif_IRhot |
---|
91 | !=================================================================================================== |
---|
92 | ! |
---|
93 | !=================================================================================================== |
---|
94 | ! function Agrif_Rhot |
---|
95 | ! |
---|
96 | !> Returns the time refinement factor of the current grid. |
---|
97 | !--------------------------------------------------------------------------------------------------- |
---|
98 | function Agrif_Rhot ( ) result( rhot ) |
---|
99 | !--------------------------------------------------------------------------------------------------- |
---|
100 | real :: rhot |
---|
101 | ! |
---|
102 | rhot = float(Agrif_IRhot()) |
---|
103 | !--------------------------------------------------------------------------------------------------- |
---|
104 | end function Agrif_Rhot |
---|
105 | !=================================================================================================== |
---|
106 | ! |
---|
107 | !=================================================================================================== |
---|
108 | ! function Agrif_Parent_IRhot |
---|
109 | ! |
---|
110 | !> Returns the time refinement factor of the parent of the current grid. |
---|
111 | !--------------------------------------------------------------------------------------------------- |
---|
112 | function Agrif_Parent_IRhot ( ) result( irhot ) |
---|
113 | !--------------------------------------------------------------------------------------------------- |
---|
114 | integer :: i, irhot |
---|
115 | ! |
---|
116 | irhot = 1 |
---|
117 | ! |
---|
118 | do i = 1,Agrif_Probdim |
---|
119 | irhot = max(irhot, Agrif_Curgrid % parent % timeref(i)) |
---|
120 | enddo |
---|
121 | !--------------------------------------------------------------------------------------------------- |
---|
122 | end function Agrif_Parent_IRhot |
---|
123 | !=================================================================================================== |
---|
124 | ! |
---|
125 | !=================================================================================================== |
---|
126 | ! function Agrif_Parent_Rhot |
---|
127 | ! |
---|
128 | !> Returns the time refinement factor of the parent of the current grid. |
---|
129 | !--------------------------------------------------------------------------------------------------- |
---|
130 | function Agrif_Parent_Rhot ( ) result( rhot ) |
---|
131 | !--------------------------------------------------------------------------------------------------- |
---|
132 | real :: rhot |
---|
133 | ! |
---|
134 | rhot = float(Agrif_Parent_IRhot()) |
---|
135 | !--------------------------------------------------------------------------------------------------- |
---|
136 | end function Agrif_Parent_Rhot |
---|
137 | !=================================================================================================== |
---|
138 | ! |
---|
139 | !=================================================================================================== |
---|
140 | ! function Agrif_Nbstepint |
---|
141 | ! |
---|
142 | !> function for the calculation of the coefficients used for the time interpolation |
---|
143 | !! (module #Agrif_Boundary). |
---|
144 | !--------------------------------------------------------------------------------------------------- |
---|
145 | function Agrif_Nbstepint ( ) |
---|
146 | !--------------------------------------------------------------------------------------------------- |
---|
147 | integer :: Agrif_nbstepint ! result |
---|
148 | ! |
---|
149 | Agrif_nbstepint = mod(Agrif_Curgrid % ngridstep, Agrif_iRhot()) |
---|
150 | !--------------------------------------------------------------------------------------------------- |
---|
151 | end function Agrif_Nbstepint |
---|
152 | !=================================================================================================== |
---|
153 | ! |
---|
154 | !=================================================================================================== |
---|
155 | ! function Agrif_Parent_Nbstepint |
---|
156 | ! |
---|
157 | !> function for the calculation of the coefficients used for the time interpolation |
---|
158 | !! (module #Agrif_Boundary). |
---|
159 | !--------------------------------------------------------------------------------------------------- |
---|
160 | function Agrif_Parent_Nbstepint ( ) |
---|
161 | !--------------------------------------------------------------------------------------------------- |
---|
162 | integer :: Agrif_Parent_Nbstepint ! result |
---|
163 | ! |
---|
164 | Agrif_Parent_Nbstepint = mod(Agrif_Curgrid % parent % ngridstep, int(Agrif_Parent_Rhot())) |
---|
165 | !--------------------------------------------------------------------------------------------------- |
---|
166 | end function Agrif_Parent_Nbstepint |
---|
167 | !=================================================================================================== |
---|
168 | ! |
---|
169 | !=================================================================================================== |
---|
170 | ! subroutine Agrif_InterpNearBorderX |
---|
171 | ! |
---|
172 | !> Allows to interpolate (in the x direction) on a near border of the current grid if this one |
---|
173 | !! has a common border with the root coarse grid. |
---|
174 | !--------------------------------------------------------------------------------------------------- |
---|
175 | subroutine Agrif_InterpNearBorderX ( ) |
---|
176 | !--------------------------------------------------------------------------------------------------- |
---|
177 | Agrif_Curgrid % NearRootBorder(1) = .FALSE. |
---|
178 | !--------------------------------------------------------------------------------------------------- |
---|
179 | end subroutine Agrif_InterpNearBorderX |
---|
180 | !=================================================================================================== |
---|
181 | ! |
---|
182 | !=================================================================================================== |
---|
183 | ! subroutine Agrif_InterpDistantBorderX |
---|
184 | ! |
---|
185 | !> Allows to interpolate (in the x direction) on a distant border of the current grid if this one |
---|
186 | !! has a common border with the root coarse grid. |
---|
187 | !--------------------------------------------------------------------------------------------------- |
---|
188 | subroutine Agrif_InterpDistantBorderX ( ) |
---|
189 | !--------------------------------------------------------------------------------------------------- |
---|
190 | Agrif_Curgrid % DistantRootBorder(1) = .FALSE. |
---|
191 | !--------------------------------------------------------------------------------------------------- |
---|
192 | end subroutine Agrif_InterpDistantBorderX |
---|
193 | !=================================================================================================== |
---|
194 | ! |
---|
195 | !=================================================================================================== |
---|
196 | ! subroutine Agrif_InterpNearBorderY |
---|
197 | ! |
---|
198 | !> Allows to interpolate (in the y direction) on a near border of the current grid if this one |
---|
199 | !! has a common border with the root coarse grid. |
---|
200 | !--------------------------------------------------------------------------------------------------- |
---|
201 | subroutine Agrif_InterpNearBorderY ( ) |
---|
202 | !--------------------------------------------------------------------------------------------------- |
---|
203 | Agrif_Curgrid % NearRootBorder(2) = .FALSE. |
---|
204 | !--------------------------------------------------------------------------------------------------- |
---|
205 | end subroutine Agrif_InterpNearBorderY |
---|
206 | !=================================================================================================== |
---|
207 | ! |
---|
208 | !=================================================================================================== |
---|
209 | ! subroutine Agrif_InterpDistantBorderY |
---|
210 | ! |
---|
211 | !> Allows to interpolate (in the y direction) on a distant border of the current grid if this one |
---|
212 | !! has a common border with the root coarse grid. |
---|
213 | !--------------------------------------------------------------------------------------------------- |
---|
214 | subroutine Agrif_InterpDistantBorderY ( ) |
---|
215 | !--------------------------------------------------------------------------------------------------- |
---|
216 | Agrif_Curgrid % DistantRootBorder(2) = .FALSE. |
---|
217 | !--------------------------------------------------------------------------------------------------- |
---|
218 | end subroutine Agrif_InterpDistantBorderY |
---|
219 | !=================================================================================================== |
---|
220 | ! |
---|
221 | !=================================================================================================== |
---|
222 | ! subroutine Agrif_InterpNearBorderZ |
---|
223 | ! |
---|
224 | !> Allows to interpolate (in the z direction) on a near border of the current grid if this one |
---|
225 | !! has a common border with the root coarse grid. |
---|
226 | !--------------------------------------------------------------------------------------------------- |
---|
227 | subroutine Agrif_InterpNearBorderZ ( ) |
---|
228 | !--------------------------------------------------------------------------------------------------- |
---|
229 | Agrif_Curgrid % NearRootBorder(3) = .FALSE. |
---|
230 | !--------------------------------------------------------------------------------------------------- |
---|
231 | end subroutine Agrif_InterpNearBorderZ |
---|
232 | !=================================================================================================== |
---|
233 | ! |
---|
234 | !=================================================================================================== |
---|
235 | ! subroutine Agrif_InterpDistantBorderZ |
---|
236 | ! |
---|
237 | !> Allows to interpolate (in the z direction) on a distant border of the current grid if this one |
---|
238 | !! has a common border with the root coarse grid. |
---|
239 | !--------------------------------------------------------------------------------------------------- |
---|
240 | subroutine Agrif_InterpDistantBorderZ() |
---|
241 | !--------------------------------------------------------------------------------------------------- |
---|
242 | Agrif_Curgrid % DistantRootBorder(3) = .FALSE. |
---|
243 | !--------------------------------------------------------------------------------------------------- |
---|
244 | end subroutine Agrif_InterpDistantBorderZ |
---|
245 | !=================================================================================================== |
---|
246 | ! |
---|
247 | !=================================================================================================== |
---|
248 | ! function Agrif_Parent_Nb_Step |
---|
249 | ! |
---|
250 | !> Returns the number of time steps of the parent of the current grid. |
---|
251 | !--------------------------------------------------------------------------------------------------- |
---|
252 | function Agrif_Parent_Nb_Step ( ) |
---|
253 | !--------------------------------------------------------------------------------------------------- |
---|
254 | integer :: Agrif_Parent_Nb_Step ! Result |
---|
255 | ! |
---|
256 | if (Agrif_Root()) then |
---|
257 | Agrif_Parent_Nb_Step = -1 |
---|
258 | else |
---|
259 | Agrif_Parent_Nb_Step = Agrif_Curgrid % parent % ngridstep |
---|
260 | endif |
---|
261 | !--------------------------------------------------------------------------------------------------- |
---|
262 | end function Agrif_Parent_Nb_Step |
---|
263 | !=================================================================================================== |
---|
264 | ! |
---|
265 | !=================================================================================================== |
---|
266 | ! function Agrif_Root |
---|
267 | ! |
---|
268 | !> Indicates if the current grid is or not the root grid. |
---|
269 | !--------------------------------------------------------------------------------------------------- |
---|
270 | function Agrif_Root ( ) |
---|
271 | !--------------------------------------------------------------------------------------------------- |
---|
272 | logical :: Agrif_Root ! Result |
---|
273 | ! |
---|
274 | Agrif_Root = (Agrif_Curgrid % fixedrank == 0) |
---|
275 | !--------------------------------------------------------------------------------------------------- |
---|
276 | end function Agrif_Root |
---|
277 | !=================================================================================================== |
---|
278 | ! |
---|
279 | !=================================================================================================== |
---|
280 | ! function Agrif_GrandMother |
---|
281 | ! |
---|
282 | !> Indicates if the current grid is or not the root grid. |
---|
283 | !--------------------------------------------------------------------------------------------------- |
---|
284 | function Agrif_GrandMother ( ) |
---|
285 | !--------------------------------------------------------------------------------------------------- |
---|
286 | logical :: Agrif_GrandMother ! Result |
---|
287 | ! |
---|
288 | Agrif_GrandMother = Agrif_Curgrid % grand_mother_grid |
---|
289 | !--------------------------------------------------------------------------------------------------- |
---|
290 | end function Agrif_GrandMother |
---|
291 | !=================================================================================================== |
---|
292 | ! |
---|
293 | !=================================================================================================== |
---|
294 | ! function Agrif_Parent_Root |
---|
295 | ! |
---|
296 | !> Indicates if the parent of the current grid is or not the root grid. |
---|
297 | !--------------------------------------------------------------------------------------------------- |
---|
298 | function Agrif_Parent_Root ( ) |
---|
299 | !--------------------------------------------------------------------------------------------------- |
---|
300 | logical :: Agrif_Parent_Root ! Result |
---|
301 | ! |
---|
302 | Agrif_Parent_Root = (Agrif_Curgrid % parent % fixedrank == 0) |
---|
303 | !--------------------------------------------------------------------------------------------------- |
---|
304 | end function Agrif_Parent_Root |
---|
305 | !=================================================================================================== |
---|
306 | ! |
---|
307 | !=================================================================================================== |
---|
308 | ! function Agrif_Fixed |
---|
309 | ! |
---|
310 | !> Returns the number of the current grid. |
---|
311 | !--------------------------------------------------------------------------------------------------- |
---|
312 | function Agrif_Fixed ( ) |
---|
313 | !--------------------------------------------------------------------------------------------------- |
---|
314 | integer :: Agrif_Fixed ! Result |
---|
315 | ! |
---|
316 | if (Agrif_Curgrid % fixed) then |
---|
317 | Agrif_Fixed = Agrif_Curgrid % fixedrank |
---|
318 | else |
---|
319 | Agrif_Fixed = -1 |
---|
320 | endif |
---|
321 | !--------------------------------------------------------------------------------------------------- |
---|
322 | end function Agrif_Fixed |
---|
323 | !=================================================================================================== |
---|
324 | ! |
---|
325 | !=================================================================================================== |
---|
326 | ! function Agrif_Parent_Fixed |
---|
327 | ! |
---|
328 | !> Returns the number of the parent of the current grid. |
---|
329 | !--------------------------------------------------------------------------------------------------- |
---|
330 | function Agrif_Parent_Fixed ( ) |
---|
331 | !--------------------------------------------------------------------------------------------------- |
---|
332 | integer :: Agrif_Parent_Fixed ! Result |
---|
333 | ! |
---|
334 | if (Agrif_Curgrid % parent % fixed) then |
---|
335 | Agrif_Parent_Fixed = Agrif_Curgrid % parent % fixedrank |
---|
336 | else |
---|
337 | Agrif_Parent_Fixed = 0 |
---|
338 | endif |
---|
339 | !--------------------------------------------------------------------------------------------------- |
---|
340 | end function Agrif_Parent_Fixed |
---|
341 | !=================================================================================================== |
---|
342 | ! |
---|
343 | !=================================================================================================== |
---|
344 | ! function Agrif_Is_Fixed |
---|
345 | ! |
---|
346 | !> Returns .TRUE. if the current grid is fixed. |
---|
347 | !--------------------------------------------------------------------------------------------------- |
---|
348 | function Agrif_Is_Fixed ( ) |
---|
349 | !--------------------------------------------------------------------------------------------------- |
---|
350 | logical :: Agrif_Is_Fixed ! Result |
---|
351 | ! |
---|
352 | Agrif_Is_Fixed = Agrif_Curgrid % fixed |
---|
353 | !--------------------------------------------------------------------------------------------------- |
---|
354 | end function Agrif_Is_Fixed |
---|
355 | !=================================================================================================== |
---|
356 | ! |
---|
357 | !=================================================================================================== |
---|
358 | ! function Agrif_Parent_Is_Fixed |
---|
359 | ! |
---|
360 | !> Returns .TRUE. if the parent of the current grid is fixed. |
---|
361 | !--------------------------------------------------------------------------------------------------- |
---|
362 | function Agrif_Parent_Is_Fixed ( ) |
---|
363 | !--------------------------------------------------------------------------------------------------- |
---|
364 | logical :: Agrif_Parent_Is_Fixed ! Result |
---|
365 | ! |
---|
366 | Agrif_Parent_Is_Fixed = Agrif_Curgrid % parent % fixed |
---|
367 | !--------------------------------------------------------------------------------------------------- |
---|
368 | end function Agrif_Parent_Is_Fixed |
---|
369 | !=================================================================================================== |
---|
370 | ! |
---|
371 | !=================================================================================================== |
---|
372 | ! function Agrif_CFixed |
---|
373 | ! |
---|
374 | !> Returns the number of the current grid. |
---|
375 | !--------------------------------------------------------------------------------------------------- |
---|
376 | function Agrif_CFixed ( ) |
---|
377 | !--------------------------------------------------------------------------------------------------- |
---|
378 | character(3) :: Agrif_CFixed ! Result |
---|
379 | ! |
---|
380 | character(3) :: cfixed |
---|
381 | integer :: fixed |
---|
382 | ! |
---|
383 | fixed = Agrif_Fixed() |
---|
384 | ! |
---|
385 | if (fixed /= -1) then |
---|
386 | ! |
---|
387 | if (fixed <= 9) then |
---|
388 | write(cfixed,'(i1)') fixed |
---|
389 | else |
---|
390 | write(cfixed,'(i2)') fixed |
---|
391 | endif |
---|
392 | ! |
---|
393 | Agrif_CFixed = cfixed |
---|
394 | |
---|
395 | if (associated(agrif_curgrid,agrif_coarsegrid)) then |
---|
396 | Agrif_CFixed = 'gm' |
---|
397 | endif |
---|
398 | ! |
---|
399 | else |
---|
400 | print*,'Call to Agrif_CFixed() on a moving grid' |
---|
401 | stop |
---|
402 | endif |
---|
403 | !--------------------------------------------------------------------------------------------------- |
---|
404 | end function Agrif_CFixed |
---|
405 | !=================================================================================================== |
---|
406 | ! |
---|
407 | !=================================================================================================== |
---|
408 | ! function Agrid_Parent_CFixed |
---|
409 | ! |
---|
410 | !> Returns the number of the parent of the current grid. |
---|
411 | !--------------------------------------------------------------------------------------------------- |
---|
412 | function Agrid_Parent_CFixed ( ) |
---|
413 | !--------------------------------------------------------------------------------------------------- |
---|
414 | character(3) :: Agrid_Parent_CFixed ! Result |
---|
415 | ! |
---|
416 | character(3) :: cfixed |
---|
417 | integer :: fixed |
---|
418 | ! |
---|
419 | fixed = Agrif_Parent_Fixed() |
---|
420 | ! |
---|
421 | if(fixed /= -1) then |
---|
422 | ! |
---|
423 | if (fixed <= 9) then |
---|
424 | write(cfixed,'(i1)')fixed |
---|
425 | else |
---|
426 | write(cfixed,'(i2)')fixed |
---|
427 | endif |
---|
428 | ! |
---|
429 | Agrid_Parent_CFixed=cfixed |
---|
430 | ! |
---|
431 | else |
---|
432 | print*,'Illegal call to Agrid_Parent_CFixed()' |
---|
433 | stop |
---|
434 | endif |
---|
435 | !--------------------------------------------------------------------------------------------------- |
---|
436 | end function Agrid_Parent_CFixed |
---|
437 | !=================================================================================================== |
---|
438 | ! |
---|
439 | !=================================================================================================== |
---|
440 | ! subroutine Agrif_ChildGrid_to_ParentGrid |
---|
441 | ! |
---|
442 | !> Make the pointer #Agrif_Curgrid point on the parent grid of the current grid. |
---|
443 | !--------------------------------------------------------------------------------------------------- |
---|
444 | subroutine Agrif_ChildGrid_to_ParentGrid ( ) |
---|
445 | !--------------------------------------------------------------------------------------------------- |
---|
446 | Agrif_Curgrid % parent % save_grid => Agrif_Curgrid |
---|
447 | call Agrif_Instance(Agrif_Curgrid%parent) |
---|
448 | !--------------------------------------------------------------------------------------------------- |
---|
449 | end subroutine Agrif_ChildGrid_to_ParentGrid |
---|
450 | !=================================================================================================== |
---|
451 | ! |
---|
452 | !=================================================================================================== |
---|
453 | ! subroutine Agrif_ParentGrid_to_ChildGrid |
---|
454 | ! |
---|
455 | !> Make the pointer #Agrif_Curgrid point on the child grid after having called the |
---|
456 | !! #Agrif_ChildGrid_to_ParentGrid subroutine. |
---|
457 | !--------------------------------------------------------------------------------------------------- |
---|
458 | subroutine Agrif_ParentGrid_to_ChildGrid ( ) |
---|
459 | !--------------------------------------------------------------------------------------------------- |
---|
460 | call Agrif_Instance(Agrif_Curgrid%save_grid) |
---|
461 | !--------------------------------------------------------------------------------------------------- |
---|
462 | end subroutine Agrif_ParentGrid_to_ChildGrid |
---|
463 | !=================================================================================================== |
---|
464 | ! |
---|
465 | !=================================================================================================== |
---|
466 | ! function Agrif_Get_Unit |
---|
467 | ! |
---|
468 | !> Returns a unit not connected to any file. |
---|
469 | !--------------------------------------------------------------------------------------------------- |
---|
470 | function Agrif_Get_Unit ( ) |
---|
471 | !--------------------------------------------------------------------------------------------------- |
---|
472 | integer :: Agrif_Get_Unit ! Result |
---|
473 | ! |
---|
474 | integer :: n |
---|
475 | logical :: op |
---|
476 | ! |
---|
477 | integer :: nunit |
---|
478 | integer :: iii, out, iiimax |
---|
479 | logical :: bexist |
---|
480 | integer,dimension(1:1000) :: forbiddenunit |
---|
481 | ! |
---|
482 | ! Load forbidden Unit if the file Agrif_forbidenUnit exist |
---|
483 | ! |
---|
484 | INQUIRE(file='Agrif_forbiddenUnit.txt', exist=bexist) |
---|
485 | ! |
---|
486 | if (.not. bexist) then |
---|
487 | ! File Agrif_forbiddenUnit.txt not found |
---|
488 | else |
---|
489 | nunit = 777 |
---|
490 | OPEN(nunit,file='Agrif_forbiddenUnit.txt', form='formatted', status="old") |
---|
491 | iii = 1 |
---|
492 | do while ( .TRUE. ) |
---|
493 | READ(nunit,*, end=99) forbiddenunit(iii) |
---|
494 | iii = iii + 1 |
---|
495 | enddo |
---|
496 | 99 continue |
---|
497 | iiimax = iii |
---|
498 | close(nunit) |
---|
499 | endif |
---|
500 | ! |
---|
501 | do n = 7,1000 |
---|
502 | ! |
---|
503 | INQUIRE(Unit=n,Opened=op) |
---|
504 | ! |
---|
505 | out = 0 |
---|
506 | if ( bexist .AND. (.NOT.op) ) then |
---|
507 | do iii = 1,iiimax |
---|
508 | if ( n == forbiddenunit(iii) ) out = 1 |
---|
509 | enddo |
---|
510 | endif |
---|
511 | ! |
---|
512 | if ( (.NOT.op) .AND. (out == 0) ) exit |
---|
513 | ! |
---|
514 | enddo |
---|
515 | ! |
---|
516 | Agrif_Get_Unit = n |
---|
517 | !--------------------------------------------------------------------------------------------------- |
---|
518 | end function Agrif_Get_Unit |
---|
519 | !=================================================================================================== |
---|
520 | ! |
---|
521 | !=================================================================================================== |
---|
522 | ! subroutine Agrif_Set_Extra_Boundary_Cells |
---|
523 | !--------------------------------------------------------------------------------------------------- |
---|
524 | subroutine Agrif_Set_Extra_Boundary_Cells ( nb_extra_cells ) |
---|
525 | !--------------------------------------------------------------------------------------------------- |
---|
526 | integer, intent(in) :: nb_extra_cells |
---|
527 | ! |
---|
528 | Agrif_Extra_Boundary_Cells = nb_extra_cells |
---|
529 | !--------------------------------------------------------------------------------------------------- |
---|
530 | end subroutine Agrif_Set_Extra_Boundary_Cells |
---|
531 | !=================================================================================================== |
---|
532 | ! |
---|
533 | !=================================================================================================== |
---|
534 | ! subroutine Agrif_Set_Efficiency |
---|
535 | !--------------------------------------------------------------------------------------------------- |
---|
536 | subroutine Agrif_Set_Efficiency ( eff ) |
---|
537 | !--------------------------------------------------------------------------------------------------- |
---|
538 | real, intent(in) :: eff |
---|
539 | ! |
---|
540 | if ( (eff < 0.) .OR. (eff > 1) ) then |
---|
541 | write(*,*) 'Error Efficiency should be between 0 and 1' |
---|
542 | stop |
---|
543 | else |
---|
544 | Agrif_Efficiency = eff |
---|
545 | endif |
---|
546 | !--------------------------------------------------------------------------------------------------- |
---|
547 | end subroutine Agrif_Set_Efficiency |
---|
548 | !=================================================================================================== |
---|
549 | ! |
---|
550 | !=================================================================================================== |
---|
551 | ! subroutine Agrif_Set_Regridding |
---|
552 | !--------------------------------------------------------------------------------------------------- |
---|
553 | subroutine Agrif_Set_Regridding ( regfreq ) |
---|
554 | !--------------------------------------------------------------------------------------------------- |
---|
555 | integer, intent(in) :: regfreq |
---|
556 | ! |
---|
557 | if (regfreq < 0) then |
---|
558 | write(*,*) 'Regridding frequency should be positive' |
---|
559 | stop |
---|
560 | else |
---|
561 | Agrif_Regridding = regfreq |
---|
562 | endif |
---|
563 | !--------------------------------------------------------------------------------------------------- |
---|
564 | end subroutine Agrif_Set_Regridding |
---|
565 | !=================================================================================================== |
---|
566 | ! |
---|
567 | !=================================================================================================== |
---|
568 | ! subroutine Agrif_Set_coeffref_x |
---|
569 | !--------------------------------------------------------------------------------------------------- |
---|
570 | subroutine Agrif_Set_coeffref_x ( coeffref ) |
---|
571 | !--------------------------------------------------------------------------------------------------- |
---|
572 | integer, intent(in) :: coeffref |
---|
573 | |
---|
574 | if (coeffref < 0) then |
---|
575 | write(*,*) 'Coefficient of raffinement should be positive' |
---|
576 | stop |
---|
577 | else |
---|
578 | Agrif_coeffref(1) = coeffref |
---|
579 | endif |
---|
580 | !--------------------------------------------------------------------------------------------------- |
---|
581 | end subroutine Agrif_Set_coeffref_x |
---|
582 | !=================================================================================================== |
---|
583 | ! |
---|
584 | !=================================================================================================== |
---|
585 | ! subroutine Agrif_Set_coeffref_y |
---|
586 | !--------------------------------------------------------------------------------------------------- |
---|
587 | subroutine Agrif_Set_coeffref_y ( coeffref ) |
---|
588 | !--------------------------------------------------------------------------------------------------- |
---|
589 | integer, intent(in) :: coeffref |
---|
590 | |
---|
591 | if (coeffref < 0) then |
---|
592 | write(*,*) 'Coefficient of raffinement should be positive' |
---|
593 | stop |
---|
594 | else |
---|
595 | Agrif_coeffref(2) = coeffref |
---|
596 | endif |
---|
597 | !--------------------------------------------------------------------------------------------------- |
---|
598 | end subroutine Agrif_Set_coeffref_y |
---|
599 | !=================================================================================================== |
---|
600 | ! |
---|
601 | !=================================================================================================== |
---|
602 | ! subroutine Agrif_Set_coeffref_z |
---|
603 | !--------------------------------------------------------------------------------------------------- |
---|
604 | subroutine Agrif_Set_coeffref_z ( coeffref ) |
---|
605 | !--------------------------------------------------------------------------------------------------- |
---|
606 | integer, intent(in) :: coeffref |
---|
607 | ! |
---|
608 | if (coeffref < 0) then |
---|
609 | write(*,*) 'Coefficient of raffinement should be positive' |
---|
610 | stop |
---|
611 | else |
---|
612 | Agrif_coeffref(3) = coeffref |
---|
613 | endif |
---|
614 | !--------------------------------------------------------------------------------------------------- |
---|
615 | end subroutine Agrif_Set_coeffref_z |
---|
616 | !=================================================================================================== |
---|
617 | ! |
---|
618 | !=================================================================================================== |
---|
619 | ! subroutine Agrif_Set_coeffreft_x |
---|
620 | !--------------------------------------------------------------------------------------------------- |
---|
621 | subroutine Agrif_Set_coeffreft_x ( coeffref ) |
---|
622 | !--------------------------------------------------------------------------------------------------- |
---|
623 | integer, intent(in) :: coeffref |
---|
624 | |
---|
625 | if (coeffref < 0) then |
---|
626 | write(*,*) 'Coefficient of time raffinement should be positive' |
---|
627 | stop |
---|
628 | else |
---|
629 | Agrif_coeffreft(1) = coeffref |
---|
630 | endif |
---|
631 | !--------------------------------------------------------------------------------------------------- |
---|
632 | end subroutine Agrif_Set_coeffreft_x |
---|
633 | !=================================================================================================== |
---|
634 | ! |
---|
635 | !=================================================================================================== |
---|
636 | ! subroutine Agrif_Set_coeffreft_y |
---|
637 | !--------------------------------------------------------------------------------------------------- |
---|
638 | subroutine Agrif_Set_coeffreft_y ( coeffref ) |
---|
639 | !--------------------------------------------------------------------------------------------------- |
---|
640 | integer, intent(in) :: coeffref |
---|
641 | ! |
---|
642 | if (coeffref < 0) then |
---|
643 | write(*,*) 'Coefficient of time raffinement should be positive' |
---|
644 | stop |
---|
645 | else |
---|
646 | Agrif_coeffreft(2) = coeffref |
---|
647 | endif |
---|
648 | !--------------------------------------------------------------------------------------------------- |
---|
649 | end subroutine Agrif_Set_coeffreft_y |
---|
650 | !=================================================================================================== |
---|
651 | ! |
---|
652 | !=================================================================================================== |
---|
653 | ! subroutine Agrif_Set_coeffreft_z |
---|
654 | !--------------------------------------------------------------------------------------------------- |
---|
655 | subroutine Agrif_Set_coeffreft_z ( coeffref ) |
---|
656 | !--------------------------------------------------------------------------------------------------- |
---|
657 | integer, intent(in) :: coeffref |
---|
658 | |
---|
659 | if (coeffref < 0) then |
---|
660 | write(*,*)'Coefficient of time raffinement should be positive' |
---|
661 | stop |
---|
662 | else |
---|
663 | Agrif_coeffreft(3) = coeffref |
---|
664 | endif |
---|
665 | !--------------------------------------------------------------------------------------------------- |
---|
666 | end subroutine Agrif_Set_coeffreft_z |
---|
667 | !=================================================================================================== |
---|
668 | ! |
---|
669 | !=================================================================================================== |
---|
670 | ! subroutine Agrif_Set_Minwidth |
---|
671 | !--------------------------------------------------------------------------------------------------- |
---|
672 | subroutine Agrif_Set_Minwidth ( coefminwidth ) |
---|
673 | !--------------------------------------------------------------------------------------------------- |
---|
674 | integer, intent(in) :: coefminwidth |
---|
675 | ! |
---|
676 | if (coefminwidth < 0) then |
---|
677 | write(*,*)'Coefficient of Minwidth should be positive' |
---|
678 | stop |
---|
679 | else |
---|
680 | Agrif_Minwidth = coefminwidth |
---|
681 | endif |
---|
682 | !--------------------------------------------------------------------------------------------------- |
---|
683 | end subroutine Agrif_Set_Minwidth |
---|
684 | !=================================================================================================== |
---|
685 | ! |
---|
686 | !=================================================================================================== |
---|
687 | ! subroutine Agrif_Set_Rafmax |
---|
688 | !--------------------------------------------------------------------------------------------------- |
---|
689 | subroutine Agrif_Set_Rafmax ( coefrafmax ) |
---|
690 | !--------------------------------------------------------------------------------------------------- |
---|
691 | integer, intent(in) :: coefrafmax |
---|
692 | ! |
---|
693 | integer :: i |
---|
694 | real :: res |
---|
695 | ! |
---|
696 | if (coefrafmax < 0) then |
---|
697 | write(*,*)'Coefficient of should be positive' |
---|
698 | stop |
---|
699 | else |
---|
700 | res = 1. |
---|
701 | do i = 1,coefrafmax-1 |
---|
702 | res = res * FLOAT(Agrif_coeffref(1)) |
---|
703 | enddo |
---|
704 | if ( res == 0 ) res = 1 |
---|
705 | Agrif_Mind(1) = 1. / res |
---|
706 | ! |
---|
707 | res = 1. |
---|
708 | do i = 1,coefrafmax-1 |
---|
709 | res = res * FLOAT(Agrif_coeffref(2)) |
---|
710 | enddo |
---|
711 | if ( res == 0 ) res = 1 |
---|
712 | Agrif_Mind(2) = 1. / res |
---|
713 | ! |
---|
714 | res = 1. |
---|
715 | do i = 1,coefrafmax-1 |
---|
716 | res = res * FLOAT(Agrif_coeffref(3)) |
---|
717 | enddo |
---|
718 | if ( res == 0 ) res = 1 |
---|
719 | Agrif_Mind(3) = 1. / res |
---|
720 | ! |
---|
721 | endif |
---|
722 | !--------------------------------------------------------------------------------------------------- |
---|
723 | end subroutine Agrif_Set_Rafmax |
---|
724 | !=================================================================================================== |
---|
725 | ! |
---|
726 | !=================================================================================================== |
---|
727 | ! subroutine Agrif_Set_MaskMaxSearch |
---|
728 | !--------------------------------------------------------------------------------------------------- |
---|
729 | subroutine Agrif_Set_MaskMaxSearch ( mymaxsearch ) |
---|
730 | !--------------------------------------------------------------------------------------------------- |
---|
731 | integer, intent(in) :: mymaxsearch |
---|
732 | ! |
---|
733 | MaxSearch = mymaxsearch |
---|
734 | !--------------------------------------------------------------------------------------------------- |
---|
735 | end subroutine Agrif_Set_MaskMaxSearch |
---|
736 | !=================================================================================================== |
---|
737 | ! |
---|
738 | !=================================================================================================== |
---|
739 | ! function Agrif_Level |
---|
740 | !--------------------------------------------------------------------------------------------------- |
---|
741 | function Agrif_Level ( ) |
---|
742 | !--------------------------------------------------------------------------------------------------- |
---|
743 | integer :: Agrif_Level ! Result |
---|
744 | ! |
---|
745 | Agrif_Level = Agrif_Curgrid % level |
---|
746 | !--------------------------------------------------------------------------------------------------- |
---|
747 | end function Agrif_Level |
---|
748 | !=================================================================================================== |
---|
749 | ! |
---|
750 | !=================================================================================================== |
---|
751 | ! function Agrif_MaxLevel |
---|
752 | !--------------------------------------------------------------------------------------------------- |
---|
753 | function Agrif_MaxLevel ( ) |
---|
754 | !--------------------------------------------------------------------------------------------------- |
---|
755 | integer :: Agrif_MaxLevel ! Result |
---|
756 | ! |
---|
757 | Agrif_MaxLevel = Agrif_MaxLevelLoc |
---|
758 | !--------------------------------------------------------------------------------------------------- |
---|
759 | end function Agrif_MaxLevel |
---|
760 | !=================================================================================================== |
---|
761 | ! |
---|
762 | !=================================================================================================== |
---|
763 | ! function Agrif_GridAllocation_is_done |
---|
764 | !--------------------------------------------------------------------------------------------------- |
---|
765 | function Agrif_GridAllocation_is_done ( ) result(isdone) |
---|
766 | !--------------------------------------------------------------------------------------------------- |
---|
767 | logical :: isdone |
---|
768 | ! |
---|
769 | isdone = Agrif_Curgrid % allocation_is_done |
---|
770 | !--------------------------------------------------------------------------------------------------- |
---|
771 | end function Agrif_GridAllocation_is_done |
---|
772 | !=================================================================================================== |
---|
773 | ! |
---|
774 | |
---|
775 | function Agrif_Parent_Real_4(real_variable) result(real_variable_parent) |
---|
776 | real(KIND=4) :: real_variable |
---|
777 | real(KIND=4) :: real_variable_parent |
---|
778 | |
---|
779 | integer :: i |
---|
780 | logical :: i_found |
---|
781 | |
---|
782 | i_found = .FALSE. |
---|
783 | |
---|
784 | do i=1,Agrif_NbVariables(2) |
---|
785 | if (LOC(real_variable) == LOC(agrif_curgrid%tabvars_r(i)%array0)) then |
---|
786 | real_variable_parent = agrif_curgrid%tabvars_r(i)%parent_var%array0 |
---|
787 | i_found = .TRUE. |
---|
788 | EXIT |
---|
789 | endif |
---|
790 | enddo |
---|
791 | |
---|
792 | IF (.NOT.i_found) THEN |
---|
793 | do i=1,Agrif_NbVariables(2) |
---|
794 | if (LOC(real_variable) == LOC(agrif_curgrid%tabvars_r(i)%sarray0)) then |
---|
795 | real_variable_parent = agrif_curgrid%tabvars_r(i)%parent_var%sarray0 |
---|
796 | i_found = .TRUE. |
---|
797 | EXIT |
---|
798 | endif |
---|
799 | enddo |
---|
800 | ENDIF |
---|
801 | |
---|
802 | if (.NOT.i_found) STOP 'Agrif_Parent_Real_4 : Variable not found' |
---|
803 | |
---|
804 | end function Agrif_Parent_Real_4 |
---|
805 | |
---|
806 | function Agrif_Parent_Real_8(real_variable) result(real_variable_parent) |
---|
807 | real(KIND=8) :: real_variable |
---|
808 | real(KIND=8) :: real_variable_parent |
---|
809 | |
---|
810 | integer :: i |
---|
811 | logical :: i_found |
---|
812 | |
---|
813 | i_found = .FALSE. |
---|
814 | |
---|
815 | do i=1,Agrif_NbVariables(2) |
---|
816 | if (LOC(real_variable) == LOC(agrif_curgrid%tabvars_r(i)%array0)) then |
---|
817 | real_variable_parent = agrif_curgrid%tabvars_r(i)%parent_var%array0 |
---|
818 | i_found = .TRUE. |
---|
819 | EXIT |
---|
820 | endif |
---|
821 | enddo |
---|
822 | |
---|
823 | IF (.NOT.i_found) THEN |
---|
824 | do i=1,Agrif_NbVariables(2) |
---|
825 | if (LOC(real_variable) == LOC(agrif_curgrid%tabvars_r(i)%darray0)) then |
---|
826 | real_variable_parent = agrif_curgrid%tabvars_r(i)%parent_var%darray0 |
---|
827 | i_found = .TRUE. |
---|
828 | EXIT |
---|
829 | endif |
---|
830 | enddo |
---|
831 | ENDIF |
---|
832 | |
---|
833 | if (.NOT.i_found) STOP 'Agrif_Parent_Real_8 : Variable not found' |
---|
834 | |
---|
835 | end function Agrif_Parent_Real_8 |
---|
836 | |
---|
837 | function Agrif_Parent_Integer(integer_variable) result(integer_variable_parent) |
---|
838 | integer :: integer_variable |
---|
839 | integer :: integer_variable_parent |
---|
840 | |
---|
841 | integer :: i |
---|
842 | logical :: i_found |
---|
843 | |
---|
844 | i_found = .FALSE. |
---|
845 | |
---|
846 | do i=1,Agrif_NbVariables(4) |
---|
847 | if (LOC(integer_variable) == LOC(agrif_curgrid%tabvars_i(i)%iarray0)) then |
---|
848 | integer_variable_parent = agrif_curgrid%tabvars_i(i)%parent_var%iarray0 |
---|
849 | i_found = .TRUE. |
---|
850 | EXIT |
---|
851 | endif |
---|
852 | enddo |
---|
853 | |
---|
854 | if (.NOT.i_found) STOP 'Agrif_Parent : Variable not found' |
---|
855 | |
---|
856 | end function Agrif_Parent_Integer |
---|
857 | |
---|
858 | function Agrif_Parent_Character(character_variable) result(character_variable_parent) |
---|
859 | character(*) :: character_variable |
---|
860 | character(len(character_variable)) :: character_variable_parent |
---|
861 | |
---|
862 | integer :: i |
---|
863 | logical :: i_found |
---|
864 | |
---|
865 | i_found = .FALSE. |
---|
866 | |
---|
867 | do i=1,Agrif_NbVariables(1) |
---|
868 | if (LOC(character_variable) == LOC(agrif_curgrid%tabvars_c(i)%carray0)) then |
---|
869 | character_variable_parent = agrif_curgrid%tabvars_c(i)%parent_var%carray0 |
---|
870 | i_found = .TRUE. |
---|
871 | EXIT |
---|
872 | endif |
---|
873 | enddo |
---|
874 | |
---|
875 | if (.NOT.i_found) STOP 'Agrif_Parent : Variable not found' |
---|
876 | |
---|
877 | end function Agrif_Parent_Character |
---|
878 | |
---|
879 | function Agrif_Parent_Logical(logical_variable) result(logical_variable_parent) |
---|
880 | logical :: logical_variable |
---|
881 | logical :: logical_variable_parent |
---|
882 | |
---|
883 | integer :: i |
---|
884 | logical :: i_found |
---|
885 | |
---|
886 | i_found = .FALSE. |
---|
887 | |
---|
888 | do i=1,Agrif_NbVariables(3) |
---|
889 | if (LOC(logical_variable) == LOC(agrif_curgrid%tabvars_l(i)%larray0)) then |
---|
890 | logical_variable_parent = agrif_curgrid%tabvars_l(i)%parent_var%larray0 |
---|
891 | i_found = .TRUE. |
---|
892 | EXIT |
---|
893 | endif |
---|
894 | enddo |
---|
895 | |
---|
896 | if (.NOT.i_found) STOP 'Agrif_Parent : Variable not found' |
---|
897 | |
---|
898 | end function Agrif_Parent_Logical |
---|
899 | |
---|
900 | function Agrif_Irhox() result(i_val) |
---|
901 | integer :: i_val |
---|
902 | i_val = agrif_curgrid%spaceref(1) |
---|
903 | end function Agrif_Irhox |
---|
904 | |
---|
905 | function Agrif_Irhoy() result(i_val) |
---|
906 | integer :: i_val |
---|
907 | i_val = agrif_curgrid%spaceref(2) |
---|
908 | end function Agrif_Irhoy |
---|
909 | |
---|
910 | function Agrif_Irhoz() result(i_val) |
---|
911 | integer :: i_val |
---|
912 | i_val = agrif_curgrid%spaceref(3) |
---|
913 | end function Agrif_Irhoz |
---|
914 | |
---|
915 | function Agrif_NearCommonBorderX() result(l_val) |
---|
916 | logical :: l_val |
---|
917 | l_val = agrif_curgrid%nearRootBorder(1) |
---|
918 | end function Agrif_NearCommonBorderX |
---|
919 | |
---|
920 | function Agrif_NearCommonBorderY() result(l_val) |
---|
921 | logical :: l_val |
---|
922 | l_val = agrif_curgrid%nearRootBorder(2) |
---|
923 | end function Agrif_NearCommonBorderY |
---|
924 | |
---|
925 | function Agrif_NearCommonBorderZ() result(l_val) |
---|
926 | logical :: l_val |
---|
927 | l_val = agrif_curgrid%nearRootBorder(3) |
---|
928 | end function Agrif_NearCommonBorderZ |
---|
929 | |
---|
930 | function Agrif_DistantCommonBorderX() result(l_val) |
---|
931 | logical :: l_val |
---|
932 | l_val = agrif_curgrid%DistantRootBorder(1) |
---|
933 | end function Agrif_DistantCommonBorderX |
---|
934 | |
---|
935 | function Agrif_DistantCommonBorderY() result(l_val) |
---|
936 | logical :: l_val |
---|
937 | l_val = agrif_curgrid%DistantRootBorder(2) |
---|
938 | end function Agrif_DistantCommonBorderY |
---|
939 | |
---|
940 | function Agrif_DistantCommonBorderZ() result(l_val) |
---|
941 | logical :: l_val |
---|
942 | l_val = agrif_curgrid%DistantRootBorder(3) |
---|
943 | end function Agrif_DistantCommonBorderZ |
---|
944 | |
---|
945 | function Agrif_Ix() result(i_val) |
---|
946 | integer :: i_val |
---|
947 | i_val = agrif_curgrid%ix(1) |
---|
948 | end function Agrif_Ix |
---|
949 | |
---|
950 | function Agrif_Iy() result(i_val) |
---|
951 | integer :: i_val |
---|
952 | i_val = agrif_curgrid%ix(2) |
---|
953 | end function Agrif_Iy |
---|
954 | |
---|
955 | function Agrif_Iz() result(i_val) |
---|
956 | integer :: i_val |
---|
957 | i_val = agrif_curgrid%ix(3) |
---|
958 | end function Agrif_Iz |
---|
959 | |
---|
960 | function Agrif_Get_grid_id() result(i_val) |
---|
961 | integer :: i_val |
---|
962 | i_val = agrif_curgrid % grid_id |
---|
963 | end function Agrif_Get_grid_id |
---|
964 | |
---|
965 | function Agrif_Get_parent_id() result(i_val) |
---|
966 | integer :: i_val |
---|
967 | i_val = agrif_curgrid % parent % grid_id |
---|
968 | end function Agrif_Get_parent_id |
---|
969 | |
---|
970 | function Agrif_rhox() result(r_val) |
---|
971 | real :: r_val |
---|
972 | r_val = real(agrif_curgrid%spaceref(1)) |
---|
973 | end function Agrif_rhox |
---|
974 | |
---|
975 | function Agrif_rhoy() result(r_val) |
---|
976 | real :: r_val |
---|
977 | r_val = real(agrif_curgrid%spaceref(2)) |
---|
978 | end function Agrif_rhoy |
---|
979 | |
---|
980 | function Agrif_rhoz() result(r_val) |
---|
981 | real :: r_val |
---|
982 | r_val = real(agrif_curgrid%spaceref(3)) |
---|
983 | end function Agrif_rhoz |
---|
984 | |
---|
985 | function Agrif_Nb_Step() result(i_val) |
---|
986 | integer :: i_val |
---|
987 | i_val = agrif_curgrid%ngridstep |
---|
988 | end function Agrif_Nb_Step |
---|
989 | |
---|
990 | function Agrif_Nb_Fine_Grids() result(i_val) |
---|
991 | integer :: i_val |
---|
992 | i_val = Agrif_nbfixedgrids |
---|
993 | end function Agrif_Nb_Fine_Grids |
---|
994 | |
---|
995 | end module Agrif_CurgridFunctions |
---|