-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathgr_initParameshArrays.F90
143 lines (131 loc) · 5.13 KB
/
gr_initParameshArrays.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
!!****if* source/Grid/GridMain/AMR/Paramesh4/PM4_package/gr_initParameshArrays
!! NOTICE
!! This file derived from PARAMESH - an adaptive mesh library.
!! Copyright (C) 2003, 2004 United States Government as represented by the
!! National Aeronautics and Space Administration, Goddard Space Flight
!! Center. All Rights Reserved.
!! Copyright 2022 UChicago Argonne, LLC and contributors
!!
!! Use of the PARAMESH software is governed by the terms of the
!! usage agreement which can be found in the file
!! 'PARAMESH_USERS_AGREEMENT' in the main paramesh directory.
!!
!! NAME
!!
!! gr_initParameshArrays
!!
!! SYNOPSIS
!!
!! call gr_initParameshArrays(logical(IN) :: restart,
!! integer(IN) :: xlboundary,
!! integer(IN) :: xrboundary,
!! integer(IN) :: ylboundary,
!! integer(IN) :: yrboundary,
!! integer(IN) :: zlboundary,
!! integer(IN) :: zrboundary
!! )
!!
!! DESCRIPTION
!!
!! Perform early initialization of some Grid data structures.
!!
!! This routine prepares the Grid for being filled with
!! meaningful data.
!!
!! ARGUMENTS
!!
!! restart - Is the grid being prepared for initialization with
!! data from a checkpoint file?
!! xlboundary - boundary condition type of outer domain boundary in lower X direction.
!! xrboundary - boundary condition type of outer domain boundary in upper X direction.
!! ylboundary - boundary condition type of outer domain boundary in lower Y direction.
!! yrboundary - boundary condition type of outer domain boundary in upper Y direction.
!! zlboundary - boundary condition type of outer domain boundary in lower Z direction.
!! zrboundary - boundary condition type of outer domain boundary in upper Z direction.
!!
!!
!!
!! HISTORY
!!
!! 2003 - 2022 Adapted for FLASH and Flash-X U of Chicago
!!***
subroutine gr_initParameshArrays(restart,&
& xlboundary, xrboundary, &
& ylboundary, yrboundary, &
& zlboundary, zrboundary)
use paramesh_dimensions
use physicaldata
use workspace
use tree
use paramesh_mpi_interfaces, ONLY : mpi_amr_global_domain_limits,&
amr_morton_process, &
mpi_amr_boundary_block_info
use paramesh_interfaces, ONLY : amr_refine_derefine
use Grid_data, ONLY : gr_meshMe, gr_meshNumProcs
use gr_specificData, ONLY : gr_gidIsValid
use gr_interface, ONLY : gr_pmIoTreeMetadataIsValid
use Logfile_interface, ONLY : Logfile_stampMessage
use Driver_interface, only: Driver_abort
use Simulation_interface, ONLY : Simulation_mapIntToStr
implicit none
#include "constants.h"
#include "Simulation.h"
logical,intent(IN) :: restart
integer,intent(IN) :: xlboundary, xrboundary
integer,intent(IN) :: ylboundary, yrboundary
integer,intent(IN) :: zlboundary, zrboundary
integer :: i
character(len=4) :: vname
#if defined(BITTREE)
call mpi_amr_global_domain_limits
call amr_build_bittree()
if(.not. restart) then
call amr_reorder_grid()
call amr_refine_derefine()
end if
#else
! Make sure that blocks produced by divide_domain are in strict
! morton order
if(restart) then
call mpi_amr_global_domain_limits()
if (.NOT. (gr_pmIoTreeMetadataIsValid() .OR. gr_gidIsValid)) then
! This will happen in particular if we are restarting from a
! checkpoint that was written by a run that used the Amrex
! Grid implementation.
! If only the surr_blks data was bad, then we can recover
! by a amr_morton_process call below and so do not need to
! abort here.
call Logfile_stampMessage("Invalid Grid tree metadata after reading checkpoint;&
& you may need to recompile with the Bittree feature if you want to restart&
& from this checkpoint with a PARAMESH Grid.")
call Driver_abort("The Grid tree metadata in the checkpoint is insufficient&
& for restarting with this Grid implementation.")
end if
else
call amr_refine_derefine()
call mpi_amr_global_domain_limits
call amr_reorder_grid()
end if
#endif
call gr_initParameshDomainBboxes(xlboundary, xrboundary, &
& ylboundary, yrboundary, &
& zlboundary, zrboundary)
call amr_morton_process()
if(restart) then
grid_analysed_mpi=1
#ifdef CALL_BOUNDARY_BLOCK_INFO
call mpi_amr_boundary_block_info(gr_meshMe,gr_meshNumProcs)
#endif
end if
! reset for quadratic interpolation
interp_mask_work(:) = 1
! AH: do not reset interp_mask_unk, but warn if not one
!interp_mask_unk(:) = 1
do i = 1, nvar
if (interp_mask_unk(i) .NE. 1) then
call Simulation_mapIntToStr(i,vname,MAPBLOCK_UNK)
call Logfile_stampMessage("WARNING: interp_mask_unk is not 1 for variable "//vname)
end if
end do
return
end subroutine gr_initParameshArrays