[4777] | 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 Agrif_Util |
---|
| 24 | !! |
---|
| 25 | !! This module contains the two procedures called in the main program : |
---|
| 26 | !! - #Agrif_Init_Grids allows the initialization of the root coarse grid |
---|
| 27 | !! - #Agrif_Step allows the creation of the grid hierarchy and the management of the time integration. |
---|
| 28 | ! |
---|
| 29 | module Agrif_Util |
---|
| 30 | ! |
---|
| 31 | use Agrif_Clustering |
---|
| 32 | use Agrif_BcFunction |
---|
| 33 | use Agrif_seq |
---|
| 34 | ! |
---|
| 35 | implicit none |
---|
| 36 | ! |
---|
| 37 | abstract interface |
---|
| 38 | subroutine step_proc() |
---|
| 39 | end subroutine step_proc |
---|
| 40 | end interface |
---|
| 41 | ! |
---|
| 42 | contains |
---|
| 43 | ! |
---|
| 44 | !=================================================================================================== |
---|
| 45 | ! subroutine Agrif_Step |
---|
| 46 | ! |
---|
| 47 | !> Creates the grid hierarchy and manages the time integration procedure. |
---|
| 48 | !> It is called in the main program. |
---|
| 49 | !> Calls subroutines #Agrif_Regrid and #Agrif_Integrate. |
---|
| 50 | !--------------------------------------------------------------------------------------------------- |
---|
| 51 | subroutine Agrif_Step ( procname ) |
---|
| 52 | !--------------------------------------------------------------------------------------------------- |
---|
| 53 | procedure(step_proc) :: procname !< subroutine to call on each grid |
---|
| 54 | type(agrif_grid), pointer :: ref_grid |
---|
| 55 | ! |
---|
| 56 | ! Set the clustering variables |
---|
| 57 | call Agrif_clustering_def() |
---|
| 58 | ! |
---|
| 59 | ! Creation and initialization of the grid hierarchy |
---|
| 60 | if ( Agrif_USE_ONLY_FIXED_GRIDS == 1 ) then |
---|
| 61 | ! |
---|
| 62 | if ( (Agrif_Mygrid % ngridstep == 0) .AND. (.not. Agrif_regrid_has_been_done) ) then |
---|
| 63 | call Agrif_Regrid() |
---|
| 64 | Agrif_regrid_has_been_done = .TRUE. |
---|
| 65 | endif |
---|
| 66 | ! |
---|
| 67 | else |
---|
| 68 | ! |
---|
| 69 | if (mod(Agrif_Mygrid % ngridstep,Agrif_Regridding) == 0) then |
---|
| 70 | call Agrif_Regrid() |
---|
| 71 | endif |
---|
| 72 | ! |
---|
| 73 | endif |
---|
| 74 | ! |
---|
| 75 | ! Time integration of the grid hierarchy |
---|
| 76 | if (agrif_coarse) then |
---|
| 77 | ref_grid => agrif_coarsegrid |
---|
| 78 | else |
---|
| 79 | ref_grid => agrif_mygrid |
---|
| 80 | endif |
---|
| 81 | if ( Agrif_Parallel_sisters ) then |
---|
| 82 | call Agrif_Integrate_Parallel(ref_grid,procname) |
---|
| 83 | else |
---|
| 84 | call Agrif_Integrate(ref_grid,procname) |
---|
| 85 | endif |
---|
| 86 | ! |
---|
| 87 | if ( ref_grid%child_list%nitems > 0 ) call Agrif_Instance(ref_grid) |
---|
| 88 | !--------------------------------------------------------------------------------------------------- |
---|
| 89 | end subroutine Agrif_Step |
---|
| 90 | !=================================================================================================== |
---|
| 91 | ! |
---|
| 92 | !=================================================================================================== |
---|
| 93 | ! subroutine Agrif_Step_Child |
---|
| 94 | ! |
---|
| 95 | !> Apply 'procname' to each grid of the hierarchy |
---|
| 96 | !--------------------------------------------------------------------------------------------------- |
---|
| 97 | subroutine Agrif_Step_Child ( procname ) |
---|
| 98 | !--------------------------------------------------------------------------------------------------- |
---|
| 99 | procedure(step_proc) :: procname !< subroutine to call on each grid |
---|
| 100 | ! |
---|
| 101 | if ( Agrif_Parallel_sisters ) then |
---|
| 102 | call Agrif_Integrate_Child_Parallel(Agrif_Mygrid,procname) |
---|
| 103 | else |
---|
| 104 | call Agrif_Integrate_Child(Agrif_Mygrid,procname) |
---|
| 105 | endif |
---|
| 106 | ! |
---|
| 107 | if ( Agrif_Mygrid%child_list%nitems > 0 ) call Agrif_Instance(Agrif_Mygrid) |
---|
| 108 | !--------------------------------------------------------------------------------------------------- |
---|
| 109 | end subroutine Agrif_Step_Child |
---|
| 110 | !=================================================================================================== |
---|
| 111 | ! |
---|
| 112 | !=================================================================================================== |
---|
| 113 | ! subroutine Agrif_Regrid |
---|
| 114 | ! |
---|
| 115 | !> Creates the grid hierarchy from fixed grids and adaptive mesh refinement. |
---|
| 116 | !--------------------------------------------------------------------------------------------------- |
---|
| 117 | subroutine Agrif_Regrid ( procname ) |
---|
| 118 | !--------------------------------------------------------------------------------------------------- |
---|
| 119 | procedure(init_proc), optional :: procname !< Initialisation subroutine (Default: Agrif_InitValues) |
---|
| 120 | ! |
---|
| 121 | type(Agrif_Rectangle), pointer :: coarsegrid_fixed |
---|
| 122 | type(Agrif_Rectangle), pointer :: coarsegrid_moving |
---|
| 123 | integer :: i, j |
---|
| 124 | integer :: nunit |
---|
| 125 | logical :: BEXIST |
---|
| 126 | TYPE(Agrif_Rectangle) :: newrect ! Pointer on a new grid |
---|
| 127 | integer :: is_coarse, rhox, rhoy, rhoz, rhot |
---|
| 128 | ! |
---|
| 129 | if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) & |
---|
| 130 | call Agrif_detect_all(Agrif_Mygrid) ! Detection of areas to be refined |
---|
| 131 | ! |
---|
| 132 | allocate(coarsegrid_fixed) |
---|
| 133 | allocate(coarsegrid_moving) |
---|
| 134 | ! |
---|
| 135 | if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) & |
---|
| 136 | call Agrif_Cluster_All(Agrif_Mygrid,coarsegrid_moving) ! Clustering |
---|
| 137 | ! |
---|
| 138 | if ( Agrif_USE_FIXED_GRIDS == 1 .OR. Agrif_USE_ONLY_FIXED_GRIDS == 1 ) then |
---|
| 139 | ! |
---|
| 140 | if (Agrif_Mygrid % ngridstep == 0) then |
---|
| 141 | ! |
---|
| 142 | nunit = Agrif_Get_Unit() |
---|
| 143 | open(nunit, file='AGRIF_FixedGrids.in', form='formatted', status="old", ERR=99) |
---|
| 144 | if (agrif_coarse) then ! SKIP the coarse grid declaration |
---|
| 145 | if (Agrif_Probdim == 3) then |
---|
| 146 | read(nunit,*) is_coarse, rhox, rhoy, rhoz, rhot |
---|
| 147 | elseif (Agrif_Probdim == 2) then |
---|
| 148 | read(nunit,*) is_coarse, rhox, rhoy, rhot |
---|
| 149 | elseif (Agrif_Probdim == 2) then |
---|
| 150 | read(nunit,*) is_coarse, rhox, rhot |
---|
| 151 | endif |
---|
| 152 | endif |
---|
| 153 | ! Creation of the grid hierarchy from the Agrif_FixedGrids.in file |
---|
| 154 | do i = 1,Agrif_Probdim |
---|
| 155 | coarsegrid_fixed % imin(i) = 1 |
---|
| 156 | coarsegrid_fixed % imax(i) = Agrif_Mygrid % nb(i) + 1 |
---|
| 157 | enddo |
---|
| 158 | j = 1 |
---|
| 159 | call Agrif_Read_Fix_Grd(coarsegrid_fixed,j,nunit) |
---|
| 160 | close(nunit) |
---|
| 161 | ! |
---|
| 162 | call Agrif_gl_clear(Agrif_oldmygrid) |
---|
| 163 | ! |
---|
| 164 | ! Creation of the grid hierarchy from coarsegrid_fixed |
---|
| 165 | call Agrif_Create_Grids(Agrif_Mygrid,coarsegrid_fixed) |
---|
| 166 | |
---|
| 167 | else |
---|
| 168 | call Agrif_gl_copy(Agrif_oldmygrid, Agrif_Mygrid % child_list) |
---|
| 169 | endif |
---|
| 170 | else |
---|
| 171 | call Agrif_gl_copy(Agrif_oldmygrid, Agrif_Mygrid % child_list) |
---|
| 172 | call Agrif_gl_clear(Agrif_Mygrid % child_list) |
---|
| 173 | endif |
---|
| 174 | ! |
---|
| 175 | if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) then |
---|
| 176 | ! |
---|
| 177 | call Agrif_Save_All(Agrif_oldmygrid) |
---|
| 178 | call Agrif_Free_before_All(Agrif_oldmygrid) |
---|
| 179 | ! |
---|
| 180 | ! Creation of the grid hierarchy from coarsegrid_moving |
---|
| 181 | call Agrif_Create_Grids(Agrif_Mygrid,coarsegrid_moving) |
---|
| 182 | ! |
---|
| 183 | endif |
---|
| 184 | ! |
---|
| 185 | ! Initialization of the grid hierarchy by copy or interpolation |
---|
| 186 | ! |
---|
| 187 | #if defined AGRIF_MPI |
---|
| 188 | if ( Agrif_Parallel_sisters ) then |
---|
| 189 | call Agrif_Init_Hierarchy_Parallel_1(Agrif_Mygrid) |
---|
| 190 | call Agrif_Init_Hierarchy_Parallel_2(Agrif_Mygrid,procname) |
---|
| 191 | else |
---|
| 192 | call Agrif_Init_Hierarchy(Agrif_Mygrid,procname) |
---|
| 193 | endif |
---|
| 194 | #else |
---|
| 195 | call Agrif_Init_Hierarchy(Agrif_Mygrid,procname) |
---|
| 196 | #endif |
---|
| 197 | ! |
---|
| 198 | if ( Agrif_USE_ONLY_FIXED_GRIDS == 0 ) call Agrif_Free_after_All(Agrif_oldmygrid) |
---|
| 199 | ! |
---|
| 200 | Agrif_regrid_has_been_done = .TRUE. |
---|
| 201 | ! |
---|
| 202 | call Agrif_Instance( Agrif_Mygrid ) |
---|
| 203 | ! |
---|
| 204 | deallocate(coarsegrid_fixed) |
---|
| 205 | deallocate(coarsegrid_moving) |
---|
| 206 | ! |
---|
| 207 | return |
---|
| 208 | ! |
---|
| 209 | ! Opening error |
---|
| 210 | ! |
---|
| 211 | 99 INQUIRE(FILE='AGRIF_FixedGrids.in',EXIST=BEXIST) |
---|
| 212 | if (.not. BEXIST) then |
---|
| 213 | print*,'ERROR : File AGRIF_FixedGrids.in not found.' |
---|
| 214 | STOP |
---|
| 215 | else |
---|
| 216 | print*,'Error opening file AGRIF_FixedGrids.in' |
---|
| 217 | STOP |
---|
| 218 | endif |
---|
| 219 | !--------------------------------------------------------------------------------------------------- |
---|
| 220 | end subroutine Agrif_Regrid |
---|
| 221 | !=================================================================================================== |
---|
| 222 | ! |
---|
| 223 | !=================================================================================================== |
---|
| 224 | ! subroutine Agrif_detect_All |
---|
| 225 | ! |
---|
| 226 | !> Detects areas to be refined. |
---|
| 227 | !--------------------------------------------------------------------------------------------------- |
---|
| 228 | recursive subroutine Agrif_detect_all ( g ) |
---|
| 229 | !--------------------------------------------------------------------------------------------------- |
---|
| 230 | TYPE(Agrif_Grid), pointer :: g !< Pointer on the current grid |
---|
| 231 | ! |
---|
| 232 | Type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure |
---|
| 233 | integer, DIMENSION(3) :: size |
---|
| 234 | integer :: i |
---|
| 235 | real :: g_eps |
---|
| 236 | ! |
---|
| 237 | parcours => g % child_list % first |
---|
| 238 | ! |
---|
| 239 | ! To be positioned on the finer grids of the grid hierarchy |
---|
| 240 | ! |
---|
| 241 | do while (associated(parcours)) |
---|
| 242 | call Agrif_detect_all(parcours % gr) |
---|
| 243 | parcours => parcours % next |
---|
| 244 | enddo |
---|
| 245 | ! |
---|
| 246 | g_eps = huge(1.) |
---|
| 247 | do i = 1,Agrif_Probdim |
---|
| 248 | g_eps = min(g_eps, g % Agrif_dx(i)) |
---|
| 249 | enddo |
---|
| 250 | ! |
---|
| 251 | g_eps = g_eps / 100. |
---|
| 252 | ! |
---|
| 253 | if ( Agrif_Probdim == 1 ) g%tabpoint1D = 0 |
---|
| 254 | if ( Agrif_Probdim == 2 ) g%tabpoint2D = 0 |
---|
| 255 | if ( Agrif_Probdim == 3 ) g%tabpoint3D = 0 |
---|
| 256 | ! |
---|
| 257 | do i = 1,Agrif_Probdim |
---|
| 258 | if ( g%Agrif_dx(i)/Agrif_coeffref(i) < (Agrif_mind(i)-g_eps) ) return |
---|
| 259 | enddo |
---|
| 260 | ! |
---|
| 261 | call Agrif_instance(g) |
---|
| 262 | ! |
---|
| 263 | ! Detection (Agrif_detect is a users routine) |
---|
| 264 | ! |
---|
| 265 | do i = 1,Agrif_Probdim |
---|
| 266 | size(i) = g % nb(i) + 1 |
---|
| 267 | enddo |
---|
| 268 | ! |
---|
| 269 | SELECT CASE (Agrif_Probdim) |
---|
| 270 | CASE (1) |
---|
| 271 | call Agrif_detect(g%tabpoint1D,size) |
---|
| 272 | CASE (2) |
---|
| 273 | call Agrif_detect(g%tabpoint2D,size) |
---|
| 274 | CASE (3) |
---|
| 275 | call Agrif_detect(g%tabpoint3D,size) |
---|
| 276 | END SELECT |
---|
| 277 | ! |
---|
| 278 | ! Addition of the areas detected on the child grids |
---|
| 279 | ! |
---|
| 280 | parcours => g % child_list % first |
---|
| 281 | ! |
---|
| 282 | do while (associated(parcours)) |
---|
| 283 | call Agrif_Add_detected_areas(g,parcours % gr) |
---|
| 284 | parcours => parcours % next |
---|
| 285 | enddo |
---|
| 286 | !--------------------------------------------------------------------------------------------------- |
---|
| 287 | end subroutine Agrif_detect_all |
---|
| 288 | !=================================================================================================== |
---|
| 289 | ! |
---|
| 290 | !=================================================================================================== |
---|
| 291 | ! subroutine Agrif_Add_detected_areas |
---|
| 292 | ! |
---|
| 293 | !> Adds on the parent grid the areas detected on its child grids |
---|
| 294 | !--------------------------------------------------------------------------------------------------- |
---|
| 295 | subroutine Agrif_Add_detected_areas ( parentgrid, childgrid ) |
---|
| 296 | !--------------------------------------------------------------------------------------------------- |
---|
| 297 | Type(Agrif_Grid), pointer :: parentgrid |
---|
| 298 | Type(Agrif_Grid), pointer :: childgrid |
---|
| 299 | ! |
---|
| 300 | integer :: i,j,k |
---|
| 301 | ! |
---|
| 302 | do i = 1,childgrid%nb(1)+1 |
---|
| 303 | if ( Agrif_Probdim == 1 ) then |
---|
| 304 | if (childgrid%tabpoint1D(i)==1) then |
---|
| 305 | parentgrid%tabpoint1D(childgrid%ix(1)+(i-1)/Agrif_Coeffref(1)) = 1 |
---|
| 306 | endif |
---|
| 307 | else |
---|
| 308 | do j=1,childgrid%nb(2)+1 |
---|
| 309 | if (Agrif_Probdim==2) then |
---|
| 310 | if (childgrid%tabpoint2D(i,j)==1) then |
---|
| 311 | parentgrid%tabpoint2D( & |
---|
| 312 | childgrid%ix(1)+(i-1)/Agrif_Coeffref(1), & |
---|
| 313 | childgrid%ix(2)+(j-1)/Agrif_Coeffref(2)) = 1 |
---|
| 314 | endif |
---|
| 315 | else |
---|
| 316 | do k=1,childgrid%nb(3)+1 |
---|
| 317 | if (childgrid%tabpoint3D(i,j,k)==1) then |
---|
| 318 | parentgrid%tabpoint3D( & |
---|
| 319 | childgrid%ix(1)+(i-1)/Agrif_Coeffref(1), & |
---|
| 320 | childgrid%ix(2)+(j-1)/Agrif_Coeffref(2), & |
---|
| 321 | childgrid%ix(3)+(k-1)/Agrif_Coeffref(3)) = 1 |
---|
| 322 | endif |
---|
| 323 | enddo |
---|
| 324 | endif |
---|
| 325 | enddo |
---|
| 326 | endif |
---|
| 327 | enddo |
---|
| 328 | !--------------------------------------------------------------------------------------------------- |
---|
| 329 | end subroutine Agrif_Add_detected_areas |
---|
| 330 | !=================================================================================================== |
---|
| 331 | ! |
---|
| 332 | !=================================================================================================== |
---|
| 333 | ! subroutine Agrif_Free_before_All |
---|
| 334 | !--------------------------------------------------------------------------------------------------- |
---|
| 335 | recursive subroutine Agrif_Free_before_All ( gridlist ) |
---|
| 336 | !--------------------------------------------------------------------------------------------------- |
---|
| 337 | Type(Agrif_Grid_List), intent(inout) :: gridlist !< Grid list |
---|
| 338 | ! |
---|
| 339 | Type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure |
---|
| 340 | ! |
---|
| 341 | parcours => gridlist % first |
---|
| 342 | ! |
---|
| 343 | do while (associated(parcours)) |
---|
| 344 | ! |
---|
| 345 | if (.not. parcours%gr%fixed) then |
---|
| 346 | call Agrif_Free_data_before(parcours%gr) |
---|
| 347 | parcours % gr % oldgrid = .TRUE. |
---|
| 348 | endif |
---|
| 349 | ! |
---|
| 350 | call Agrif_Free_before_all (parcours % gr % child_list) |
---|
| 351 | ! |
---|
| 352 | parcours => parcours % next |
---|
| 353 | ! |
---|
| 354 | enddo |
---|
| 355 | !--------------------------------------------------------------------------------------------------- |
---|
| 356 | end subroutine Agrif_Free_before_All |
---|
| 357 | !=================================================================================================== |
---|
| 358 | ! |
---|
| 359 | !=================================================================================================== |
---|
| 360 | ! subroutine Agrif_Save_All |
---|
| 361 | !--------------------------------------------------------------------------------------------------- |
---|
| 362 | recursive subroutine Agrif_Save_All ( gridlist ) |
---|
| 363 | !--------------------------------------------------------------------------------------------------- |
---|
| 364 | type(Agrif_Grid_List), intent(inout) :: gridlist !< Grid list |
---|
| 365 | ! |
---|
| 366 | type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure |
---|
| 367 | ! |
---|
| 368 | parcours => gridlist % first |
---|
| 369 | ! |
---|
| 370 | do while (associated(parcours)) |
---|
| 371 | ! |
---|
| 372 | if (.not. parcours%gr%fixed) then |
---|
| 373 | call Agrif_Instance(parcours%gr) |
---|
| 374 | call Agrif_Before_Regridding() |
---|
| 375 | parcours % gr % oldgrid = .TRUE. |
---|
| 376 | endif |
---|
| 377 | ! |
---|
| 378 | call Agrif_Save_All(parcours % gr % child_list) |
---|
| 379 | ! |
---|
| 380 | parcours => parcours % next |
---|
| 381 | ! |
---|
| 382 | enddo |
---|
| 383 | !--------------------------------------------------------------------------------------------------- |
---|
| 384 | end subroutine Agrif_Save_All |
---|
| 385 | !=================================================================================================== |
---|
| 386 | ! |
---|
| 387 | !=================================================================================================== |
---|
| 388 | ! subroutine Agrif_Free_after_All |
---|
| 389 | !--------------------------------------------------------------------------------------------------- |
---|
| 390 | recursive subroutine Agrif_Free_after_All ( gridlist ) |
---|
| 391 | !--------------------------------------------------------------------------------------------------- |
---|
| 392 | Type(Agrif_Grid_List), intent(inout) :: gridlist !< Grid list to free |
---|
| 393 | ! |
---|
| 394 | Type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive proced |
---|
| 395 | Type(Agrif_PGrid), pointer :: preparcours |
---|
| 396 | Type(Agrif_PGrid), pointer :: preparcoursini |
---|
| 397 | ! |
---|
| 398 | allocate(preparcours) |
---|
| 399 | ! |
---|
| 400 | preparcoursini => preparcours |
---|
| 401 | ! |
---|
| 402 | nullify(preparcours % gr) |
---|
| 403 | ! |
---|
| 404 | preparcours % next => gridlist % first |
---|
| 405 | parcours => gridlist % first |
---|
| 406 | ! |
---|
| 407 | do while (associated(parcours)) |
---|
| 408 | ! |
---|
| 409 | if ( (.NOT. parcours%gr % fixed) .AND. (parcours%gr % oldgrid) ) then |
---|
| 410 | call Agrif_Free_data_after(parcours%gr) |
---|
| 411 | endif |
---|
| 412 | ! |
---|
| 413 | call Agrif_Free_after_all( parcours%gr % child_list ) |
---|
| 414 | ! |
---|
| 415 | if (parcours % gr % oldgrid) then |
---|
| 416 | deallocate(parcours % gr) |
---|
| 417 | preparcours % next => parcours % next |
---|
| 418 | deallocate(parcours) |
---|
| 419 | parcours => preparcours % next |
---|
| 420 | else |
---|
| 421 | preparcours => preparcours % next |
---|
| 422 | parcours => parcours % next |
---|
| 423 | endif |
---|
| 424 | ! |
---|
| 425 | enddo |
---|
| 426 | ! |
---|
| 427 | deallocate(preparcoursini) |
---|
| 428 | !--------------------------------------------------------------------------------------------------- |
---|
| 429 | end subroutine Agrif_Free_after_All |
---|
| 430 | !=================================================================================================== |
---|
| 431 | ! |
---|
| 432 | !=================================================================================================== |
---|
| 433 | ! subroutine Agrif_Integrate |
---|
| 434 | ! |
---|
| 435 | !> Manages the time integration of the grid hierarchy. |
---|
| 436 | !! Recursive subroutine and call on subroutines Agrif_Init::Agrif_Instance and #Agrif_Step |
---|
| 437 | !--------------------------------------------------------------------------------------------------- |
---|
| 438 | recursive subroutine Agrif_Integrate ( g, procname ) |
---|
| 439 | !--------------------------------------------------------------------------------------------------- |
---|
| 440 | type(Agrif_Grid), pointer :: g !< Pointer on the current grid |
---|
| 441 | procedure(step_proc) :: procname !< Subroutine to call on each grid |
---|
| 442 | ! |
---|
| 443 | type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure |
---|
| 444 | integer :: nbt ! Number of time steps of the current grid |
---|
| 445 | integer :: i, k |
---|
| 446 | ! |
---|
| 447 | ! Instanciation of the variables of the current grid |
---|
| 448 | ! if ( g % fixedrank /= 0 ) then |
---|
| 449 | call Agrif_Instance(g) |
---|
| 450 | ! endif |
---|
| 451 | ! |
---|
| 452 | ! One step on the current grid |
---|
| 453 | ! |
---|
| 454 | call procname () |
---|
| 455 | ! |
---|
| 456 | ! Number of time steps on the current grid |
---|
| 457 | ! |
---|
| 458 | g%ngridstep = g % ngridstep + 1 |
---|
| 459 | parcours => g % child_list % first |
---|
| 460 | ! |
---|
| 461 | ! Recursive procedure for the time integration of the grid hierarchy |
---|
| 462 | do while (associated(parcours)) |
---|
| 463 | ! |
---|
| 464 | ! Instanciation of the variables of the current grid |
---|
| 465 | call Agrif_Instance(parcours % gr) |
---|
| 466 | ! |
---|
| 467 | ! Number of time steps |
---|
| 468 | nbt = 1 |
---|
| 469 | do i = 1,Agrif_Probdim |
---|
| 470 | nbt = max(nbt, parcours % gr % timeref(i)) |
---|
| 471 | enddo |
---|
| 472 | ! |
---|
| 473 | do k = 1,nbt |
---|
| 474 | call Agrif_Integrate(parcours % gr, procname) |
---|
| 475 | enddo |
---|
| 476 | ! |
---|
| 477 | parcours => parcours % next |
---|
| 478 | ! |
---|
| 479 | enddo |
---|
| 480 | !--------------------------------------------------------------------------------------------------- |
---|
| 481 | end subroutine Agrif_Integrate |
---|
| 482 | !=================================================================================================== |
---|
| 483 | ! |
---|
| 484 | !=================================================================================================== |
---|
| 485 | ! subroutine Agrif_Integrate_Parallel |
---|
| 486 | ! |
---|
| 487 | !> Manages the time integration of the grid hierarchy in parallel |
---|
| 488 | !! Recursive subroutine and call on subroutines Agrif_Init::Agrif_Instance and #Agrif_Step |
---|
| 489 | !--------------------------------------------------------------------------------------------------- |
---|
| 490 | recursive subroutine Agrif_Integrate_Parallel ( g, procname ) |
---|
| 491 | !--------------------------------------------------------------------------------------------------- |
---|
| 492 | type(Agrif_Grid), pointer :: g !< Pointer on the current grid |
---|
| 493 | procedure(step_proc) :: procname !< Subroutine to call on each grid |
---|
| 494 | ! |
---|
| 495 | #if defined AGRIF_MPI |
---|
| 496 | type(Agrif_PGrid), pointer :: gridp ! Pointer for the recursive procedure |
---|
| 497 | integer :: nbt ! Number of time steps of the current grid |
---|
| 498 | integer :: i, k, is |
---|
| 499 | ! |
---|
| 500 | ! Instanciation of the variables of the current grid |
---|
| 501 | if ( g % fixedrank /= 0 ) then |
---|
| 502 | call Agrif_Instance(g) |
---|
| 503 | endif |
---|
| 504 | ! |
---|
| 505 | ! One step on the current grid |
---|
| 506 | call procname () |
---|
| 507 | ! |
---|
| 508 | ! Number of time steps on the current grid |
---|
| 509 | g % ngridstep = g % ngridstep + 1 |
---|
| 510 | ! |
---|
| 511 | ! Continue only if the grid has defined sequences of child integrations. |
---|
| 512 | if ( .not. associated(g % child_seq) ) return |
---|
| 513 | ! |
---|
| 514 | do is = 1, g % child_seq % nb_seqs |
---|
| 515 | ! |
---|
| 516 | ! For each sequence, a given processor does integrate only on grid. |
---|
| 517 | gridp => Agrif_seq_select_child(g,is) |
---|
| 518 | ! |
---|
| 519 | ! Instanciation of the variables of the current grid |
---|
| 520 | call Agrif_Instance(gridp % gr) |
---|
| 521 | ! |
---|
| 522 | ! Number of time steps |
---|
| 523 | nbt = 1 |
---|
| 524 | do i = 1,Agrif_Probdim |
---|
| 525 | nbt = max(nbt, gridp % gr % timeref(i)) |
---|
| 526 | enddo |
---|
| 527 | ! |
---|
| 528 | do k = 1,nbt |
---|
| 529 | call Agrif_Integrate_Parallel(gridp % gr, procname) |
---|
| 530 | enddo |
---|
| 531 | ! |
---|
| 532 | enddo |
---|
| 533 | #else |
---|
| 534 | call Agrif_Integrate( g, procname ) |
---|
| 535 | #endif |
---|
| 536 | !--------------------------------------------------------------------------------------------------- |
---|
| 537 | end subroutine Agrif_Integrate_Parallel |
---|
| 538 | !=================================================================================================== |
---|
| 539 | ! |
---|
| 540 | ! |
---|
| 541 | !=================================================================================================== |
---|
| 542 | ! subroutine Agrif_Integrate_ChildGrids |
---|
| 543 | ! |
---|
| 544 | !> Manages the time integration of the grid hierarchy. |
---|
| 545 | !! Call the subroutine procname on each child grid of the current grid |
---|
| 546 | !--------------------------------------------------------------------------------------------------- |
---|
| 547 | recursive subroutine Agrif_Integrate_ChildGrids ( procname ) |
---|
| 548 | !--------------------------------------------------------------------------------------------------- |
---|
| 549 | procedure(step_proc) :: procname !< Subroutine to call on each grid |
---|
| 550 | ! |
---|
| 551 | type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure |
---|
| 552 | integer :: nbt ! Number of time steps of the current grid |
---|
| 553 | integer :: i, k, is |
---|
| 554 | type(Agrif_Grid) , pointer :: save_grid |
---|
| 555 | type(Agrif_PGrid), pointer :: gridp ! Pointer for the recursive procedure |
---|
| 556 | |
---|
| 557 | save_grid => Agrif_Curgrid |
---|
| 558 | |
---|
| 559 | ! Number of time steps on the current grid |
---|
| 560 | save_grid % ngridstep = save_grid % ngridstep + 1 |
---|
| 561 | |
---|
| 562 | #ifdef AGRIF_MPI |
---|
| 563 | if ( .not. Agrif_Parallel_sisters ) then |
---|
| 564 | #endif |
---|
| 565 | parcours => save_grid % child_list % first |
---|
| 566 | ! |
---|
| 567 | ! Recursive procedure for the time integration of the grid hierarchy |
---|
| 568 | do while (associated(parcours)) |
---|
| 569 | ! |
---|
| 570 | ! Instanciation of the variables of the current grid |
---|
| 571 | call Agrif_Instance(parcours % gr) |
---|
| 572 | ! |
---|
| 573 | ! Number of time steps |
---|
| 574 | nbt = 1 |
---|
| 575 | do i = 1,Agrif_Probdim |
---|
| 576 | nbt = max(nbt, parcours % gr % timeref(i)) |
---|
| 577 | enddo |
---|
| 578 | ! |
---|
| 579 | do k = 1,nbt |
---|
| 580 | call procname() |
---|
| 581 | enddo |
---|
| 582 | ! |
---|
| 583 | parcours => parcours % next |
---|
| 584 | ! |
---|
| 585 | enddo |
---|
| 586 | |
---|
| 587 | #ifdef AGRIF_MPI |
---|
| 588 | else |
---|
| 589 | #endif |
---|
| 590 | ! Continue only if the grid has defined sequences of child integrations. |
---|
| 591 | if ( .not. associated(save_grid % child_seq) ) return |
---|
| 592 | ! |
---|
| 593 | do is = 1, save_grid % child_seq % nb_seqs |
---|
| 594 | ! |
---|
| 595 | ! For each sequence, a given processor does integrate only on grid. |
---|
| 596 | gridp => Agrif_seq_select_child(save_grid,is) |
---|
| 597 | ! |
---|
| 598 | ! Instanciation of the variables of the current grid |
---|
| 599 | call Agrif_Instance(gridp % gr) |
---|
| 600 | ! |
---|
| 601 | ! Number of time steps |
---|
| 602 | nbt = 1 |
---|
| 603 | do i = 1,Agrif_Probdim |
---|
| 604 | nbt = max(nbt, gridp % gr % timeref(i)) |
---|
| 605 | enddo |
---|
| 606 | ! |
---|
| 607 | do k = 1,nbt |
---|
| 608 | call procname() |
---|
| 609 | enddo |
---|
| 610 | ! |
---|
| 611 | enddo |
---|
| 612 | #ifdef AGRIF_MPI |
---|
| 613 | endif |
---|
| 614 | #endif |
---|
| 615 | |
---|
| 616 | call Agrif_Instance(save_grid) |
---|
| 617 | |
---|
| 618 | !--------------------------------------------------------------------------------------------------- |
---|
| 619 | end subroutine Agrif_Integrate_ChildGrids |
---|
| 620 | !=================================================================================================== |
---|
| 621 | !=================================================================================================== |
---|
| 622 | ! subroutine Agrif_Integrate_Child |
---|
| 623 | ! |
---|
| 624 | !> Manages the time integration of the grid hierarchy. |
---|
| 625 | !! Recursive subroutine and call on subroutines Agrif_Instance & Agrif_Step. |
---|
| 626 | !--------------------------------------------------------------------------------------------------- |
---|
| 627 | recursive subroutine Agrif_Integrate_Child ( g, procname ) |
---|
| 628 | !--------------------------------------------------------------------------------------------------- |
---|
| 629 | type(Agrif_Grid), pointer :: g !< Pointer on the current grid |
---|
| 630 | procedure(step_proc) :: procname !< Subroutine to call on each grid |
---|
| 631 | ! |
---|
| 632 | type(Agrif_PGrid), pointer :: parcours ! Pointer for the recursive procedure |
---|
| 633 | ! |
---|
| 634 | ! One step on the current grid |
---|
| 635 | ! |
---|
| 636 | call procname () |
---|
| 637 | ! |
---|
| 638 | ! Number of time steps on the current grid |
---|
| 639 | ! |
---|
| 640 | parcours => g % child_list % first |
---|
| 641 | ! |
---|
| 642 | ! Recursive procedure for the time integration of the grid hierarchy |
---|
| 643 | do while (associated(parcours)) |
---|
| 644 | ! |
---|
| 645 | ! Instanciation of the variables of the current grid |
---|
| 646 | call Agrif_Instance(parcours % gr) |
---|
| 647 | call Agrif_Integrate_Child (parcours % gr, procname) |
---|
| 648 | parcours => parcours % next |
---|
| 649 | ! |
---|
| 650 | enddo |
---|
| 651 | !--------------------------------------------------------------------------------------------------- |
---|
| 652 | end subroutine Agrif_Integrate_Child |
---|
| 653 | !=================================================================================================== |
---|
| 654 | ! |
---|
| 655 | !=================================================================================================== |
---|
| 656 | ! subroutine Agrif_Integrate_Child_Parallel |
---|
| 657 | ! |
---|
| 658 | !> Manages the time integration of the grid hierarchy. |
---|
| 659 | !! Recursive subroutine and call on subroutines Agrif_Instance & Agrif_Step. |
---|
| 660 | !--------------------------------------------------------------------------------------------------- |
---|
| 661 | recursive subroutine Agrif_Integrate_Child_Parallel ( g, procname ) |
---|
| 662 | !--------------------------------------------------------------------------------------------------- |
---|
| 663 | type(Agrif_Grid), pointer :: g !< Pointer on the current grid |
---|
| 664 | procedure(step_proc) :: procname !< Subroutine to call on each grid |
---|
| 665 | ! |
---|
| 666 | #if defined AGRIF_MPI |
---|
| 667 | type(Agrif_PGrid), pointer :: gridp ! Pointer for the recursive procedure |
---|
| 668 | integer :: is |
---|
| 669 | ! |
---|
| 670 | ! Instanciation of the variables of the current grid |
---|
| 671 | call Agrif_Instance(g) |
---|
| 672 | ! |
---|
| 673 | ! One step on the current grid |
---|
| 674 | call procname () |
---|
| 675 | ! |
---|
| 676 | ! Continue only if the grid has defined sequences of child integrations. |
---|
| 677 | if ( .not. associated(g % child_seq) ) return |
---|
| 678 | ! |
---|
| 679 | do is = 1, g % child_seq % nb_seqs |
---|
| 680 | ! |
---|
| 681 | ! For each sequence, a given processor does integrate only on grid. |
---|
| 682 | gridp => Agrif_seq_select_child(g,is) |
---|
| 683 | call Agrif_Integrate_Child_Parallel(gridp % gr, procname) |
---|
| 684 | ! |
---|
| 685 | enddo |
---|
| 686 | ! |
---|
| 687 | call Agrif_Instance(g) |
---|
| 688 | #else |
---|
| 689 | call Agrif_Integrate_Child( g, procname ) |
---|
| 690 | #endif |
---|
| 691 | !--------------------------------------------------------------------------------------------------- |
---|
| 692 | end subroutine Agrif_Integrate_Child_Parallel |
---|
| 693 | !=================================================================================================== |
---|
| 694 | ! |
---|
| 695 | !=================================================================================================== |
---|
| 696 | ! subroutine Agrif_Init_Grids |
---|
| 697 | ! |
---|
| 698 | !> Initializes the root coarse grid pointed by Agrif_Mygrid. It is called in the main program. |
---|
| 699 | !--------------------------------------------------------------------------------------------------- |
---|
| 700 | subroutine Agrif_Init_Grids ( procname1, procname2 ) |
---|
| 701 | !--------------------------------------------------------------------------------------------------- |
---|
| 702 | procedure(typedef_proc), optional :: procname1 !< (Default: Agrif_probdim_modtype_def) |
---|
| 703 | procedure(alloc_proc), optional :: procname2 !< (Default: Agrif_Allocationcalls) |
---|
| 704 | ! |
---|
| 705 | integer :: i, ierr_allocate, nunit |
---|
| 706 | integer :: is_coarse, rhox,rhoy,rhoz,rhot |
---|
| 707 | logical :: BEXIST |
---|
| 708 | ! |
---|
| 709 | if (present(procname1)) Then |
---|
| 710 | call procname1() |
---|
| 711 | else |
---|
| 712 | call Agrif_probdim_modtype_def() |
---|
| 713 | endif |
---|
| 714 | ! |
---|
| 715 | |
---|
| 716 | ! TEST FOR COARSE GRID (GRAND MOTHER GRID) in AGRIF_FixedGrids.in |
---|
| 717 | nunit = Agrif_Get_Unit() |
---|
| 718 | open(nunit, file='AGRIF_FixedGrids.in', form='formatted', status="old", ERR=98) |
---|
| 719 | if (Agrif_Probdim == 3) then |
---|
| 720 | read(nunit,*) is_coarse, rhox, rhoy, rhoz, rhot |
---|
| 721 | elseif (Agrif_Probdim == 2) then |
---|
| 722 | read(nunit,*) is_coarse, rhox, rhoy, rhot |
---|
| 723 | elseif (Agrif_Probdim == 2) then |
---|
| 724 | read(nunit,*) is_coarse, rhox, rhot |
---|
| 725 | endif |
---|
| 726 | if (is_coarse == -1) then |
---|
| 727 | agrif_coarse = .TRUE. |
---|
| 728 | if (Agrif_Probdim == 3) then |
---|
| 729 | coarse_spaceref(1:3)=(/rhox,rhoy,rhoz/) |
---|
| 730 | elseif (Agrif_Probdim == 2) then |
---|
| 731 | coarse_spaceref(1:2)=(/rhox,rhoy/) |
---|
| 732 | elseif (Agrif_Probdim == 2) then |
---|
| 733 | coarse_spaceref(1:1)=(/rhox/) |
---|
| 734 | endif |
---|
| 735 | coarse_timeref(1:Agrif_Probdim) = rhot |
---|
| 736 | endif |
---|
| 737 | close(nunit) |
---|
| 738 | |
---|
| 739 | Agrif_UseSpecialValue = .FALSE. |
---|
| 740 | Agrif_UseSpecialValueFineGrid = .FALSE. |
---|
| 741 | Agrif_SpecialValue = 0. |
---|
| 742 | Agrif_SpecialValueFineGrid = 0. |
---|
| 743 | ! |
---|
| 744 | allocate(Agrif_Mygrid) |
---|
| 745 | allocate(Agrif_OldMygrid) |
---|
| 746 | ! |
---|
| 747 | ! Space and time refinement factors are set to 1 on the root grid |
---|
| 748 | ! |
---|
| 749 | do i = 1,Agrif_Probdim |
---|
| 750 | Agrif_Mygrid % spaceref(i) = coarse_spaceref(i) |
---|
| 751 | Agrif_Mygrid % timeref(i) = coarse_timeref(i) |
---|
| 752 | enddo |
---|
| 753 | ! |
---|
| 754 | ! Initialization of the number of time steps |
---|
| 755 | Agrif_Mygrid % ngridstep = 0 |
---|
| 756 | Agrif_Mygrid % grid_id = 0 |
---|
| 757 | ! |
---|
| 758 | ! No parent grid for the root coarse grid |
---|
| 759 | nullify(Agrif_Mygrid % parent) |
---|
| 760 | ! |
---|
| 761 | ! Initialization of the minimum positions, global abscissa and space steps |
---|
| 762 | do i = 1, Agrif_Probdim |
---|
| 763 | Agrif_Mygrid % ix(i) = 1 |
---|
| 764 | Agrif_Mygrid % Agrif_x(i) = 0. |
---|
| 765 | Agrif_Mygrid % Agrif_dx(i) = 1./Agrif_Mygrid % spaceref(i) |
---|
| 766 | Agrif_Mygrid % Agrif_dt(i) = 1./Agrif_Mygrid % timeref(i) |
---|
| 767 | ! Borders of the root coarse grid |
---|
| 768 | Agrif_Mygrid % NearRootBorder(i) = .true. |
---|
| 769 | Agrif_Mygrid % DistantRootBorder(i) = .true. |
---|
| 770 | enddo |
---|
| 771 | ! |
---|
| 772 | ! The root coarse grid is a fixed grid |
---|
| 773 | Agrif_Mygrid % fixed = .TRUE. |
---|
| 774 | ! Level of the root grid |
---|
| 775 | Agrif_Mygrid % level = 0 |
---|
| 776 | ! Maximum level in the hierarchy |
---|
| 777 | Agrif_MaxLevelLoc = 0 |
---|
| 778 | ! |
---|
| 779 | ! Number of the grid pointed by Agrif_Mygrid (root coarse grid) |
---|
| 780 | Agrif_Mygrid % rank = 1 |
---|
| 781 | ! |
---|
| 782 | ! Number of the root grid as a fixed grid |
---|
| 783 | Agrif_Mygrid % fixedrank = 0 |
---|
| 784 | ! |
---|
| 785 | ! Initialization of some fields of the root grid variables |
---|
| 786 | ierr_allocate = 0 |
---|
| 787 | if( Agrif_NbVariables(0) > 0 ) allocate(Agrif_Mygrid % tabvars(Agrif_NbVariables(0)),stat=ierr_allocate) |
---|
| 788 | if( Agrif_NbVariables(1) > 0 ) allocate(Agrif_Mygrid % tabvars_c(Agrif_NbVariables(1)),stat=ierr_allocate) |
---|
| 789 | if( Agrif_NbVariables(2) > 0 ) allocate(Agrif_Mygrid % tabvars_r(Agrif_NbVariables(2)),stat=ierr_allocate) |
---|
| 790 | if( Agrif_NbVariables(3) > 0 ) allocate(Agrif_Mygrid % tabvars_l(Agrif_NbVariables(3)),stat=ierr_allocate) |
---|
| 791 | if( Agrif_NbVariables(4) > 0 ) allocate(Agrif_Mygrid % tabvars_i(Agrif_NbVariables(4)),stat=ierr_allocate) |
---|
| 792 | if (ierr_allocate /= 0) THEN |
---|
| 793 | STOP "*** ERROR WHEN ALLOCATING TABVARS ***" |
---|
| 794 | endif |
---|
| 795 | ! |
---|
| 796 | ! Initialization of the other fields of the root grid variables (number of |
---|
| 797 | ! cells, positions, number and type of its dimensions, ...) |
---|
| 798 | call Agrif_Instance(Agrif_Mygrid) |
---|
| 799 | call Agrif_Set_numberofcells(Agrif_Mygrid) |
---|
| 800 | ! |
---|
| 801 | ! Allocation of the array containing the values of the grid variables |
---|
| 802 | call Agrif_Allocation(Agrif_Mygrid, procname2) |
---|
| 803 | call Agrif_initialisations(Agrif_Mygrid) |
---|
| 804 | ! |
---|
| 805 | ! Total number of fixed grids |
---|
| 806 | Agrif_nbfixedgrids = 0 |
---|
| 807 | |
---|
| 808 | ! If a grand mother grid is declared |
---|
| 809 | |
---|
| 810 | if (agrif_coarse) then |
---|
| 811 | allocate(Agrif_Coarsegrid) |
---|
| 812 | |
---|
| 813 | Agrif_Coarsegrid % ngridstep = 0 |
---|
| 814 | Agrif_Coarsegrid % grid_id = -9999 |
---|
| 815 | |
---|
| 816 | do i = 1, Agrif_Probdim |
---|
| 817 | Agrif_Coarsegrid%spaceref(i) = coarse_spaceref(i) |
---|
| 818 | Agrif_Coarsegrid%timeref(i) = coarse_timeref(i) |
---|
| 819 | Agrif_Coarsegrid % ix(i) = 1 |
---|
| 820 | Agrif_Coarsegrid % Agrif_x(i) = 0. |
---|
| 821 | Agrif_Coarsegrid % Agrif_dx(i) = 1. |
---|
| 822 | Agrif_Coarsegrid % Agrif_dt(i) = 1. |
---|
| 823 | ! Borders of the root coarse grid |
---|
| 824 | Agrif_Coarsegrid % NearRootBorder(i) = .true. |
---|
| 825 | Agrif_Coarsegrid % DistantRootBorder(i) = .true. |
---|
| 826 | Agrif_Coarsegrid % nb(i) =Agrif_mygrid%nb(i) / coarse_spaceref(i) |
---|
| 827 | enddo |
---|
| 828 | |
---|
| 829 | ! The root coarse grid is a fixed grid |
---|
| 830 | Agrif_Coarsegrid % fixed = .TRUE. |
---|
| 831 | ! Level of the root grid |
---|
| 832 | Agrif_Coarsegrid % level = -1 |
---|
| 833 | |
---|
| 834 | Agrif_Coarsegrid % grand_mother_grid = .true. |
---|
| 835 | |
---|
| 836 | ! Number of the grid pointed by Agrif_Mygrid (root coarse grid) |
---|
| 837 | Agrif_Coarsegrid % rank = -9999 |
---|
| 838 | ! |
---|
| 839 | ! Number of the root grid as a fixed grid |
---|
| 840 | Agrif_Coarsegrid % fixedrank = -9999 |
---|
| 841 | |
---|
| 842 | Agrif_Mygrid%parent => Agrif_Coarsegrid |
---|
| 843 | |
---|
| 844 | ! Not used but required to prevent seg fault |
---|
| 845 | Agrif_Coarsegrid%parent => Agrif_Mygrid |
---|
| 846 | |
---|
| 847 | call Agrif_Create_Var(Agrif_Coarsegrid) |
---|
| 848 | |
---|
| 849 | ! Reset to null |
---|
| 850 | Nullify(Agrif_Coarsegrid%parent) |
---|
| 851 | |
---|
| 852 | Agrif_Coarsegrid%child_list%nitems = 1 |
---|
| 853 | allocate(Agrif_Coarsegrid%child_list%first) |
---|
| 854 | allocate(Agrif_Coarsegrid%child_list%last) |
---|
| 855 | Agrif_Coarsegrid%child_list%first%gr => Agrif_Mygrid |
---|
| 856 | Agrif_Coarsegrid%child_list%last%gr => Agrif_Mygrid |
---|
| 857 | |
---|
| 858 | endif |
---|
| 859 | |
---|
| 860 | return |
---|
| 861 | |
---|
| 862 | 98 INQUIRE(FILE='AGRIF_FixedGrids.in',EXIST=BEXIST) |
---|
| 863 | if (.not. BEXIST) then |
---|
| 864 | print*,'ERROR : File AGRIF_FixedGrids.in not found.' |
---|
| 865 | STOP |
---|
| 866 | else |
---|
| 867 | print*,'Error opening file AGRIF_FixedGrids.in' |
---|
| 868 | STOP |
---|
| 869 | endif |
---|
| 870 | |
---|
| 871 | !--------------------------------------------------------------------------------------------------- |
---|
| 872 | end subroutine Agrif_Init_Grids |
---|
| 873 | !=================================================================================================== |
---|
| 874 | ! |
---|
| 875 | !=================================================================================================== |
---|
| 876 | ! subroutine Agrif_Deallocation |
---|
| 877 | ! |
---|
| 878 | !> Deallocates all data arrays. |
---|
| 879 | !--------------------------------------------------------------------------------------------------- |
---|
| 880 | subroutine Agrif_Deallocation |
---|
| 881 | !--------------------------------------------------------------------------------------------------- |
---|
| 882 | integer :: nb |
---|
| 883 | type(Agrif_Variable), pointer :: var |
---|
| 884 | type(Agrif_Variable_c), pointer :: var_c |
---|
| 885 | type(Agrif_Variable_l), pointer :: var_l |
---|
| 886 | type(Agrif_Variable_i), pointer :: var_i |
---|
| 887 | ! |
---|
| 888 | do nb = 1,Agrif_NbVariables(0) |
---|
| 889 | ! |
---|
| 890 | var => Agrif_Mygrid % tabvars(nb) |
---|
| 891 | ! |
---|
| 892 | if ( allocated(var % array1) ) deallocate(var % array1) |
---|
| 893 | if ( allocated(var % array2) ) deallocate(var % array2) |
---|
| 894 | if ( allocated(var % array3) ) deallocate(var % array3) |
---|
| 895 | if ( allocated(var % array4) ) deallocate(var % array4) |
---|
| 896 | if ( allocated(var % array5) ) deallocate(var % array5) |
---|
| 897 | if ( allocated(var % array6) ) deallocate(var % array6) |
---|
| 898 | ! |
---|
| 899 | if ( allocated(var % sarray1) ) deallocate(var % sarray1) |
---|
| 900 | if ( allocated(var % sarray2) ) deallocate(var % sarray2) |
---|
| 901 | if ( allocated(var % sarray3) ) deallocate(var % sarray3) |
---|
| 902 | if ( allocated(var % sarray4) ) deallocate(var % sarray4) |
---|
| 903 | if ( allocated(var % sarray5) ) deallocate(var % sarray5) |
---|
| 904 | if ( allocated(var % sarray6) ) deallocate(var % sarray6) |
---|
| 905 | ! |
---|
| 906 | if ( allocated(var % darray1) ) deallocate(var % darray1) |
---|
| 907 | if ( allocated(var % darray2) ) deallocate(var % darray2) |
---|
| 908 | if ( allocated(var % darray3) ) deallocate(var % darray3) |
---|
| 909 | if ( allocated(var % darray4) ) deallocate(var % darray4) |
---|
| 910 | if ( allocated(var % darray5) ) deallocate(var % darray5) |
---|
| 911 | if ( allocated(var % darray6) ) deallocate(var % darray6) |
---|
| 912 | ! |
---|
| 913 | enddo |
---|
| 914 | ! |
---|
| 915 | do nb = 1,Agrif_NbVariables(1) |
---|
| 916 | ! |
---|
| 917 | var_c => Agrif_Mygrid % tabvars_c(nb) |
---|
| 918 | ! |
---|
| 919 | if ( allocated(var_c % carray1) ) deallocate(var_c % carray1) |
---|
| 920 | if ( allocated(var_c % carray2) ) deallocate(var_c % carray2) |
---|
| 921 | ! |
---|
| 922 | enddo |
---|
| 923 | |
---|
| 924 | do nb = 1,Agrif_NbVariables(3) |
---|
| 925 | ! |
---|
| 926 | var_l => Agrif_Mygrid % tabvars_l(nb) |
---|
| 927 | ! |
---|
| 928 | if ( allocated(var_l % larray1) ) deallocate(var_l % larray1) |
---|
| 929 | if ( allocated(var_l % larray2) ) deallocate(var_l % larray2) |
---|
| 930 | if ( allocated(var_l % larray3) ) deallocate(var_l % larray3) |
---|
| 931 | if ( allocated(var_l % larray4) ) deallocate(var_l % larray4) |
---|
| 932 | if ( allocated(var_l % larray5) ) deallocate(var_l % larray5) |
---|
| 933 | if ( allocated(var_l % larray6) ) deallocate(var_l % larray6) |
---|
| 934 | ! |
---|
| 935 | enddo |
---|
| 936 | ! |
---|
| 937 | do nb = 1,Agrif_NbVariables(4) |
---|
| 938 | ! |
---|
| 939 | var_i => Agrif_Mygrid % tabvars_i(nb) |
---|
| 940 | ! |
---|
| 941 | if ( allocated(var_i % iarray1) ) deallocate(var_i % iarray1) |
---|
| 942 | if ( allocated(var_i % iarray2) ) deallocate(var_i % iarray2) |
---|
| 943 | if ( allocated(var_i % iarray3) ) deallocate(var_i % iarray3) |
---|
| 944 | if ( allocated(var_i % iarray4) ) deallocate(var_i % iarray4) |
---|
| 945 | if ( allocated(var_i % iarray5) ) deallocate(var_i % iarray5) |
---|
| 946 | if ( allocated(var_i % iarray6) ) deallocate(var_i % iarray6) |
---|
| 947 | ! |
---|
| 948 | enddo |
---|
| 949 | ! |
---|
| 950 | if ( allocated(Agrif_Mygrid % tabvars) ) deallocate(Agrif_Mygrid % tabvars) |
---|
| 951 | if ( allocated(Agrif_Mygrid % tabvars_c) ) deallocate(Agrif_Mygrid % tabvars_c) |
---|
| 952 | if ( allocated(Agrif_Mygrid % tabvars_r) ) deallocate(Agrif_Mygrid % tabvars_r) |
---|
| 953 | if ( allocated(Agrif_Mygrid % tabvars_l) ) deallocate(Agrif_Mygrid % tabvars_l) |
---|
| 954 | if ( allocated(Agrif_Mygrid % tabvars_i) ) deallocate(Agrif_Mygrid % tabvars_i) |
---|
| 955 | deallocate(Agrif_Mygrid) |
---|
| 956 | !--------------------------------------------------------------------------------------------------- |
---|
| 957 | end subroutine Agrif_Deallocation |
---|
| 958 | !=================================================================================================== |
---|
| 959 | ! |
---|
| 960 | !=================================================================================================== |
---|
| 961 | ! subroutine Agrif_Step_adj |
---|
| 962 | ! |
---|
| 963 | !> creates the grid hierarchy and manages the backward time integration procedure. |
---|
| 964 | !> It is called in the main program. |
---|
| 965 | !> calls subroutines #Agrif_Regrid and #Agrif_Integrate_adj. |
---|
| 966 | !--------------------------------------------------------------------------------------------------- |
---|
| 967 | subroutine Agrif_Step_adj ( procname ) |
---|
| 968 | !--------------------------------------------------------------------------------------------------- |
---|
| 969 | procedure(step_proc) :: procname !< Subroutine to call on each grid |
---|
| 970 | ! |
---|
| 971 | ! Creation and initialization of the grid hierarchy |
---|
| 972 | ! |
---|
| 973 | ! Set the clustering variables |
---|
| 974 | call Agrif_clustering_def() |
---|
| 975 | ! |
---|
| 976 | if ( Agrif_USE_ONLY_FIXED_GRIDS .EQ. 1 ) then |
---|
| 977 | ! |
---|
| 978 | if (Agrif_Mygrid % ngridstep == 0) then |
---|
| 979 | if (.not.Agrif_regrid_has_been_done ) then |
---|
| 980 | call Agrif_Regrid() |
---|
| 981 | endif |
---|
| 982 | call Agrif_Instance(Agrif_Mygrid) |
---|
| 983 | endif |
---|
| 984 | ! |
---|
| 985 | else |
---|
| 986 | ! |
---|
| 987 | if (mod(Agrif_Mygrid % ngridstep, Agrif_Regridding) == 0) then |
---|
| 988 | call Agrif_Regrid() |
---|
| 989 | call Agrif_Instance(Agrif_Mygrid) |
---|
| 990 | endif |
---|
| 991 | ! |
---|
| 992 | endif |
---|
| 993 | ! |
---|
| 994 | ! Time integration of the grid hierarchy |
---|
| 995 | ! |
---|
| 996 | call Agrif_Integrate_adj (Agrif_Mygrid,procname) |
---|
| 997 | ! |
---|
| 998 | if ( Agrif_Mygrid % child_list % nitems > 0 ) call Agrif_Instance(Agrif_Mygrid) |
---|
| 999 | ! |
---|
| 1000 | !--------------------------------------------------------------------------------------------------- |
---|
| 1001 | end subroutine Agrif_Step_adj |
---|
| 1002 | !=================================================================================================== |
---|
| 1003 | ! |
---|
| 1004 | !=================================================================================================== |
---|
| 1005 | ! subroutine Agrif_Integrate_adj |
---|
| 1006 | ! |
---|
| 1007 | !> Manages the backward time integration of the grid hierarchy. |
---|
| 1008 | !! Recursive subroutine and call on subroutines Agrif_Init::Agrif_Instance and #Agrif_Step_adj |
---|
| 1009 | !--------------------------------------------------------------------------------------------------- |
---|
| 1010 | recursive subroutine Agrif_Integrate_adj ( g, procname ) |
---|
| 1011 | !--------------------------------------------------------------------------------------------------- |
---|
| 1012 | type(Agrif_Grid), pointer :: g !< Pointer on the current grid |
---|
| 1013 | procedure(step_proc) :: procname !< Subroutine to call on each grid |
---|
| 1014 | ! |
---|
| 1015 | type(Agrif_pgrid), pointer :: parcours ! pointer for the recursive procedure |
---|
| 1016 | integer :: nbt ! Number of time steps of the current grid |
---|
| 1017 | integer :: k |
---|
| 1018 | ! |
---|
| 1019 | ! Instanciation of the variables of the current grid |
---|
| 1020 | if ( g%fixedrank /= 0 ) then |
---|
| 1021 | call Agrif_Instance(g) |
---|
| 1022 | endif |
---|
| 1023 | ! |
---|
| 1024 | ! Number of time steps on the current grid |
---|
| 1025 | ! |
---|
| 1026 | g%ngridstep = g % ngridstep + 1 |
---|
| 1027 | parcours => g % child_list % first |
---|
| 1028 | ! |
---|
| 1029 | ! Recursive procedure for the time integration of the grid hierarchy |
---|
| 1030 | do while (associated(parcours)) |
---|
| 1031 | ! |
---|
| 1032 | ! Instanciation of the variables of the current grid |
---|
| 1033 | call Agrif_Instance(parcours % gr) |
---|
| 1034 | ! |
---|
| 1035 | ! Number of time steps |
---|
| 1036 | nbt = 1 |
---|
| 1037 | do k = 1,Agrif_Probdim |
---|
| 1038 | nbt = max(nbt, parcours % gr % timeref(k)) |
---|
| 1039 | enddo |
---|
| 1040 | ! |
---|
| 1041 | do k = nbt,1,-1 |
---|
| 1042 | call Agrif_Integrate_adj(parcours % gr, procname) |
---|
| 1043 | enddo |
---|
| 1044 | ! |
---|
| 1045 | parcours => parcours % next |
---|
| 1046 | ! |
---|
| 1047 | enddo |
---|
| 1048 | ! |
---|
| 1049 | if ( g % child_list % nitems > 0 ) call Agrif_Instance(g) |
---|
| 1050 | ! |
---|
| 1051 | ! One step on the current grid |
---|
| 1052 | call procname () |
---|
| 1053 | ! |
---|
| 1054 | end subroutine Agrif_Integrate_adj |
---|
| 1055 | !=================================================================================================== |
---|
| 1056 | ! |
---|
| 1057 | !=================================================================================================== |
---|
| 1058 | ! subroutine Agrif_Step_Child_adj |
---|
| 1059 | ! |
---|
| 1060 | !> Apply 'procname' to each grid of the hierarchy from Child to Parent |
---|
| 1061 | !--------------------------------------------------------------------------------------------------- |
---|
| 1062 | subroutine Agrif_Step_Child_adj ( procname ) |
---|
| 1063 | !--------------------------------------------------------------------------------------------------- |
---|
| 1064 | procedure(step_proc) :: procname !< Subroutine to call on each grid |
---|
| 1065 | ! |
---|
| 1066 | call Agrif_Integrate_Child_adj(Agrif_Mygrid,procname) |
---|
| 1067 | ! |
---|
| 1068 | if ( Agrif_Mygrid % child_list % nitems > 0 ) call Agrif_Instance(Agrif_Mygrid) |
---|
| 1069 | ! |
---|
| 1070 | end subroutine Agrif_Step_Child_adj |
---|
| 1071 | !=================================================================================================== |
---|
| 1072 | ! |
---|
| 1073 | !=================================================================================================== |
---|
| 1074 | ! subroutine Agrif_Integrate_Child_adj |
---|
| 1075 | ! |
---|
| 1076 | !> Manages the backward time integration of the grid hierarchy. |
---|
| 1077 | !! Recursive subroutine and call on subroutines Agrif_Init::Agrif_Instance & Agrif_Step_adj. |
---|
| 1078 | !--------------------------------------------------------------------------------------------------- |
---|
| 1079 | recursive subroutine Agrif_Integrate_Child_adj ( g, procname ) |
---|
| 1080 | !--------------------------------------------------------------------------------------------------- |
---|
| 1081 | type(Agrif_Grid),pointer :: g !< Pointer on the current grid |
---|
| 1082 | procedure(step_proc) :: procname !< Subroutine to call on each grid |
---|
| 1083 | ! |
---|
| 1084 | type(Agrif_PGrid),pointer :: parcours !< Pointer for the recursive procedure |
---|
| 1085 | ! |
---|
| 1086 | parcours => g % child_list % first |
---|
| 1087 | ! |
---|
| 1088 | ! Recursive procedure for the time integration of the grid hierarchy |
---|
| 1089 | do while (associated(parcours)) |
---|
| 1090 | ! |
---|
| 1091 | ! Instanciation of the variables of the current grid |
---|
| 1092 | call Agrif_Instance(parcours % gr) |
---|
| 1093 | call Agrif_Integrate_Child_adj(parcours % gr, procname) |
---|
| 1094 | ! |
---|
| 1095 | parcours => parcours % next |
---|
| 1096 | ! |
---|
| 1097 | enddo |
---|
| 1098 | if ( g % child_list % nitems > 0 ) call Agrif_Instance(g) |
---|
| 1099 | ! |
---|
| 1100 | ! One step on the current grid |
---|
| 1101 | call procname() |
---|
| 1102 | !--------------------------------------------------------------------------------------------------- |
---|
| 1103 | end subroutine Agrif_Integrate_Child_adj |
---|
| 1104 | !=================================================================================================== |
---|
| 1105 | ! |
---|
| 1106 | !=================================================================================================== |
---|
| 1107 | |
---|
| 1108 | end module Agrif_Util |
---|