1 |
jmc |
1.1 |
C $Header: $ |
2 |
|
|
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "CPP_OPTIONS.h" |
5 |
|
|
|
6 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
7 |
|
|
|
8 |
|
|
PROGRAM NEST_DRIVER |
9 |
|
|
C !DESCRIPTION: |
10 |
|
|
C Routine that manages the MPI communication between the CHILD |
11 |
|
|
C and PARENT models. It performs also the necessary |
12 |
|
|
C interpolations from PARENT2CHILD and CHILD2PARENT. |
13 |
|
|
C |
14 |
|
|
C ver 1.0 by G. Sannino, V. Ruggiero, A. Carillo, P. Heimbach |
15 |
|
|
C |
16 |
|
|
C First application described in: |
17 |
|
|
C Sannino G.,Herrmann, Carillo, Rupolo, Ruggiero, Artale, Heimbach, 2009: |
18 |
|
|
C An eddy-permitting model of the Mediterranean Sea with a two-way grid |
19 |
|
|
C refinement at Gibraltar. |
20 |
|
|
C Ocean Modelling, 30(1), 56-72, doi: 10.1016/j.ocemod.2009.06.002 |
21 |
|
|
C |
22 |
|
|
|
23 |
|
|
C !LOCAL INPUT VARIABLES: |
24 |
|
|
C --------------------------------------------------------------------------------- |
25 |
|
|
C NST_LEV_TOT :: Total nesting levels (1 for only one nesting) |
26 |
|
|
C NST_LEV :: Number of the actual nesting |
27 |
|
|
C NCPUs_CHLD :: Number of CPUs used for the CHILD model |
28 |
|
|
C at NST_LEV nesting level |
29 |
|
|
C NCPUs_PRNT :: Number of CPUs used for the PARENT model |
30 |
|
|
C at NST_LEV nesting level |
31 |
|
|
C nSxC,nSyC :: Domain decomposition used for CHILD |
32 |
|
|
C nSxP,nSyP :: Domain decomposition used for PARENT |
33 |
|
|
C OLX,OLY :: Domain dec. overlapping (same for both models) |
34 |
|
|
C NrP,NyP,NxP :: Dimension PARENT model |
35 |
|
|
C NrC,NyC,NxC :: Dimension CHILD model |
36 |
|
|
C WesterB :: Western (i) PARENT index where start the refinement |
37 |
|
|
C EasterB :: Eastern (i) PARENT index where finish the refinement |
38 |
|
|
C dirNEST :: Directory where are stored the geometry data of both models |
39 |
|
|
C n3dC :: number of 3-D fields sent from CHILD |
40 |
|
|
C --------------------------------------------------------------------------------- |
41 |
|
|
CEOP |
42 |
|
|
IMPLICIT NONE |
43 |
|
|
C-------------------------------------------------------- |
44 |
|
|
C INPUT VARIABLE DEFINITION |
45 |
|
|
C-------------------------------------------------------- |
46 |
|
|
INTEGER :: NST_LEV_TOT, NST_LEV, NCPUs_PRNT |
47 |
|
|
INTEGER :: Count_Lev |
48 |
|
|
PARAMETER (NST_LEV_TOT = 1) !Number of Total Nesting Levels |
49 |
|
|
PARAMETER (NST_LEV = 1) !Which level am I? |
50 |
|
|
|
51 |
|
|
INTEGER :: NCPUs_CHLD(NST_LEV_TOT) |
52 |
|
|
INTEGER :: MSTR_DRV(NST_LEV_TOT) |
53 |
|
|
INTEGER :: MSTR_PRNT(NST_LEV_TOT) |
54 |
|
|
INTEGER :: MSTR_CHLD(NST_LEV_TOT) |
55 |
|
|
|
56 |
|
|
PARAMETER (NCPUs_PRNT = 40) |
57 |
|
|
|
58 |
|
|
DATA NCPUs_CHLD / 20 / |
59 |
|
|
C-------------------------------------------------------- |
60 |
|
|
INTEGER :: nSxC,nSyC !Domain decomposition CHILD |
61 |
|
|
INTEGER :: nSxP,nSyP !Domain decomposition PARENT |
62 |
|
|
PARAMETER (nSxC = 4 , nSyC = 5) |
63 |
|
|
PARAMETER (nSxP = 8 , nSyP = 5) |
64 |
|
|
C-------------------------------------------------------- |
65 |
|
|
INTEGER :: OLY,OLX !Domain decomposition overlapping |
66 |
|
|
C !(the same for both models) |
67 |
|
|
PARAMETER (OLX = 3, OLY = 3) |
68 |
|
|
C-------------------------------------------------------- |
69 |
|
|
INTEGER :: NrP,NxP,NyP |
70 |
|
|
INTEGER :: NrC,NxC,NyC |
71 |
|
|
INTEGER :: IM_C,JM_C |
72 |
|
|
INTEGER :: IM_P,JM_P |
73 |
|
|
INTEGER :: IndI,IndJ |
74 |
|
|
INTEGER :: IndI_P(nSxP*nSyP),IndJ_P(nSxP*nSyP) |
75 |
|
|
INTEGER :: IndI_C(nSxC*nSyC),IndJ_C(nSxC*nSyC) |
76 |
|
|
|
77 |
|
|
INTEGER :: WesternB,EasternB |
78 |
|
|
C-------------------------------------------------------- |
79 |
|
|
PARAMETER (NrP=42, NyP=120,NxP = 400) !PARENT model |
80 |
|
|
PARAMETER (NrC=42, NyC=105,NxC = 140) !CHILD model |
81 |
|
|
C-------------------------------------------------------- |
82 |
|
|
PARAMETER (WesternB = 43,EasternB=90) |
83 |
|
|
C-------------------------------------------------------- |
84 |
|
|
CHARACTER :: dirNEST*80 |
85 |
|
|
C-------------------------------------------------------- |
86 |
|
|
PARAMETER (dirNEST ="/home/sannino/NESTING/") |
87 |
|
|
C-------------------------------------------------------- |
88 |
|
|
INCLUDE 'mpif.h' |
89 |
|
|
INTEGER :: ierr,rank,size,npd |
90 |
|
|
INTEGER :: irank,isize |
91 |
|
|
INTEGER :: color |
92 |
|
|
INTEGER :: istatus,NEST_comm |
93 |
|
|
INTEGER :: from,whm,status(MPI_STATUS_SIZE),st_count |
94 |
|
|
INTEGER :: I,J,K,II,JJ,Irec,III,JJJ,KK,ICONT |
95 |
|
|
INTEGER :: iVar,Indx,Jndx |
96 |
|
|
INTEGER :: J1,J2,JJ1,JJ2 |
97 |
|
|
INTEGER :: I_START,I_END,I_STEP |
98 |
|
|
REAL*4 :: XF,YF,XP1,YP1,XP2,YP2,YP3 |
99 |
|
|
REAL*8 :: TRANSPORT_WEST,TRANSPORT_EAST |
100 |
|
|
CHARACTER*10 :: c2i(30) |
101 |
|
|
C---------------------------------------------------- |
102 |
|
|
C Define PARENT Model Geometry |
103 |
|
|
C---------------------------------------------------- |
104 |
|
|
c REAL*4 Xu_P(NxP,NyP) |
105 |
|
|
REAL*4 :: Yu_P(NxP,NyP) |
106 |
|
|
REAL*4 :: Xv_P(NxP,NyP) |
107 |
|
|
REAL*4 :: Yv_P(NxP,NyP) |
108 |
|
|
REAL*4 :: Xo_P(NxP,NyP) |
109 |
|
|
REAL*4 :: Yo_P(NxP,NyP) |
110 |
|
|
REAL*4 :: Xg_P(NxP,NyP) |
111 |
|
|
REAL*4 :: Yg_P(NxP,NyP) |
112 |
|
|
REAL*4 :: hFacW_P(NxP,NyP,NrP) |
113 |
|
|
REAL*4 :: hFacS_P(NxP,NyP,NrP) |
114 |
|
|
REAL*4 :: RAC_P(NxP,NyP) |
115 |
|
|
REAL*4 :: RAW_P(NxP,NyP) |
116 |
|
|
REAL*4 :: RAS_P(NxP,NyP) |
117 |
|
|
REAL*4 :: hFacC_P(NxP,NyP,NrP) |
118 |
|
|
REAL*4 :: DEEP_P(NxP,NyP,NrP) |
119 |
|
|
REAL*4 :: INV_VOL_C_P(NxP,NyP,NrP) |
120 |
|
|
REAL*4 :: INV_VOL_S_P(NxP,NyP,NrP) |
121 |
|
|
REAL*4 :: INV_VOL_W_P(NxP,NyP,NrP) |
122 |
|
|
C---------------------------------------------------- |
123 |
|
|
C Define CHILD Model Geometry |
124 |
|
|
C---------------------------------------------------- |
125 |
|
|
REAL*4 :: Xu_C(NxC,NyC) |
126 |
|
|
REAL*4 :: Yu_C(NxC,NyC) |
127 |
|
|
REAL*4 :: Xv_C(NxC,NyC) |
128 |
|
|
REAL*4 :: Yv_C(NxC,NyC) |
129 |
|
|
REAL*4 :: Xo_C(NxC,NyC) |
130 |
|
|
REAL*4 :: Yo_C(NxC,NyC) |
131 |
|
|
REAL*4 :: Xg_C(NxC,NyC) |
132 |
|
|
REAL*4 :: Yg_C(NxC,NyC) |
133 |
|
|
REAL*4 :: hFacW_C(NxC,NyC,NrC) |
134 |
|
|
REAL*4 :: hFacS_C(NxC,NyC,NrC) |
135 |
|
|
REAL*4 :: RAC_C(NxC,NyC) |
136 |
|
|
REAL*4 :: RAW_C(NxC,NyC) |
137 |
|
|
REAL*4 :: RAS_C(NxC,NyC) |
138 |
|
|
REAL*4 :: hFacC_C(NxC,NyC,NrC) |
139 |
|
|
REAL*4 :: DEEP_C(NxC,NyC,NrC) |
140 |
|
|
C---------------------------------------------------- |
141 |
|
|
C Define relative (PARENT-->CHILD) indicies |
142 |
|
|
C---------------------------------------------------- |
143 |
|
|
INTEGER :: P2C_U(NyC) |
144 |
|
|
|
145 |
|
|
INTEGER :: P2C_linU(NyC) !Linear interp. |
146 |
|
|
INTEGER :: WO3_linU(NyC) !Linear interp. !Which Of 3 |
147 |
|
|
|
148 |
|
|
INTEGER :: P2C_linV(NyC) !Linear interp. |
149 |
|
|
INTEGER :: WO3_linV(NyC) !Linear interp. !Which Of 3 |
150 |
|
|
|
151 |
|
|
INTEGER :: P2C_V(NyC) !Linear interp. |
152 |
|
|
INTEGER :: P2C_o(NyC) !Linear interp. |
153 |
|
|
INTEGER :: P2C1_V(NyC) !BiLinear interp. |
154 |
|
|
INTEGER :: P2C2_V(NyC) !BiLinear interp. |
155 |
|
|
INTEGER :: P2C1_o(NyC) !BiLinear interp. |
156 |
|
|
INTEGER :: P2C2_o(NyC) !BiLinear interp. |
157 |
|
|
C---------------------------------------------------- |
158 |
|
|
C Define relative (CHILD-->PARENT) indicies |
159 |
|
|
C---------------------------------------------------- |
160 |
|
|
INTEGER I_C2P(9,NxP,NyP) |
161 |
|
|
INTEGER J_C2P(9,NxP,NyP) |
162 |
|
|
C---------------------------------------------------- |
163 |
|
|
C Define CHILD model variable |
164 |
|
|
C---------------------------------------------------- |
165 |
|
|
C _____________ (1) WesternB (2) EasternB |
166 |
|
|
C | |
167 |
|
|
REAL*8 :: U_C1(NyC,NrC,2) |
168 |
|
|
REAL*8 :: V_C1(NyC,NrC,2) |
169 |
|
|
REAL*8 :: T_C1(NyC,NrC,2) |
170 |
|
|
REAL*8 :: S_C1(NyC,NrC,2) |
171 |
|
|
REAL*8 :: ETA_C1(NyC,NrC,2) |
172 |
|
|
INTEGER :: MSIZE |
173 |
|
|
|
174 |
|
|
REAL*8 :: U_C2(NyC,NrC,2) |
175 |
|
|
REAL*8 :: V_C2(NyC,NrC,2) |
176 |
|
|
REAL*8 :: T_C2(NyC,NrC,2) |
177 |
|
|
REAL*8 :: S_C2(NyC,NrC,2) |
178 |
|
|
REAL*8 :: ETA_C2(NyC,NrC,2) |
179 |
|
|
|
180 |
|
|
REAL*8,allocatable :: VAR_C1(:,:,:,:) |
181 |
|
|
|
182 |
|
|
REAL*8 :: DIFF_U(NyC,NrC,2) |
183 |
|
|
REAL*8 :: DIFF_V(NyC,NrC,2) |
184 |
|
|
REAL*8 :: DIFF_T(NyC,NrC,2) |
185 |
|
|
REAL*8 :: DIFF_S(NyC,NrC,2) |
186 |
|
|
REAL*8 :: DIFF_ETA(NyC,NrC,2) |
187 |
|
|
C---------------------------------------------------- |
188 |
|
|
C Define PARENT model variable |
189 |
|
|
C---------------------------------------------------- |
190 |
|
|
REAL*8 :: VAR3D_P(NxP,NyP,NrP,15) |
191 |
|
|
REAL*8 :: VAR2D_P(NxP,NyP,4) |
192 |
|
|
|
193 |
|
|
REAL*8,allocatable :: localP3D_a(:,:,:), localP2D_a(:,:) |
194 |
|
|
|
195 |
|
|
INTEGER :: ONOFF |
196 |
|
|
INTEGER :: index_var3D,index_var2D |
197 |
|
|
C---------------------------------------------------------------| |
198 |
|
|
C (1) U || (2) V || (3) T || (4) S | |
199 |
|
|
C---------------------------------------------------------------| |
200 |
|
|
C (5) gU || (6) gV || (7) gT || (8) gS | |
201 |
|
|
C---------------------------------------------------------------| |
202 |
|
|
C (9) gUNm1 || (10) gVNm1 || (11) gTNm1 || (12) gSNm1 | |
203 |
|
|
C---------------------------------------------------------------| |
204 |
|
|
C (13) totPhiHyd || (14) IVDConvCount || | |
205 |
|
|
C---------------------------------------------------------------| |
206 |
|
|
C (15) etaN || (16) etaH || (17) phiHydLow | |
207 |
|
|
C---------------------------------------------------------------| |
208 |
|
|
C (18) etaNm1 || (19) etaHm1|| | |
209 |
|
|
C---------------------------------------------------------------| |
210 |
|
|
C---------------------------------------------------------------| |
211 |
|
|
C Define Global Variables to Exchange | |
212 |
|
|
C---------------------------------------------------------------| |
213 |
|
|
REAL*8,allocatable :: globalPA (:,:,:,:) !(6,NyP,NrP,5) |
214 |
|
|
REAL*8 :: globalP1(6,NyP,NrP) |
215 |
|
|
REAL*8 :: globalP2(6,NyP,NrP) |
216 |
|
|
REAL*8 :: globalP3(6,NyP,NrP) |
217 |
|
|
REAL*8 :: globalP4(6,NyP,NrP) |
218 |
|
|
REAL*8 :: globalP5(6,NyP,NrP) |
219 |
|
|
REAL*8 :: globalP6(6,NyP,NrP) |
220 |
|
|
REAL*8 :: globalP7(6,NyP,NrP) |
221 |
|
|
REAL*8 :: globalP8(6,NyP,NrP) |
222 |
|
|
REAL*8 :: globalP9(6,NyP,NrP) |
223 |
|
|
REAL*8 :: globalP10(6,NyP,NrP) |
224 |
|
|
REAL*8 :: globalP11(6,NyP,NrP) |
225 |
|
|
REAL*8 :: globalP12(6,NyP,NrP) |
226 |
|
|
REAL*8 :: globalP13(6,NyP,NrP) |
227 |
|
|
REAL*8 :: globalP14(6,NyP,NrP) |
228 |
|
|
|
229 |
|
|
INTEGER :: index |
230 |
|
|
C---------------------------------------------------- |
231 |
|
|
C Define Global Variables to Exchange |
232 |
|
|
C---------------------------------------------------- |
233 |
|
|
INTEGER :: n3dC |
234 |
|
|
PARAMETER ( n3dC = 15 ) |
235 |
|
|
REAL*8 :: globalC3D(NxC,NyC,NrC,n3dC) |
236 |
|
|
C |___________ 15 fields |
237 |
|
|
|
238 |
|
|
REAL*8 globalC2D(NxC,NyC,4) |
239 |
|
|
C |___________ 4 fields |
240 |
|
|
|
241 |
|
|
REAL*8,allocatable :: globalC3D_a(:,:,:,:),globalC2D_a(:,:,:) |
242 |
|
|
|
243 |
|
|
INTEGER :: indexF,index2F,index3F |
244 |
|
|
REAL*4 :: AREA_VOL |
245 |
|
|
INTEGER :: vstart,vstop,VCONT,VCONTP(0:3) |
246 |
|
|
|
247 |
|
|
C- log-file IO-unit and name: STDlog.xxxx |
248 |
|
|
INTEGER iUnit |
249 |
|
|
PARAMETER ( iUnit = 35 ) |
250 |
|
|
CHARACTER*11 fNam |
251 |
|
|
INTEGER mLoop |
252 |
|
|
C---------------------------------------------------- |
253 |
|
|
C MPI starts here |
254 |
|
|
C---------------------------------------------------- |
255 |
|
|
CALL MPI_Init(ierr) |
256 |
|
|
CALL MPI_Comm_size(MPI_COMM_WORLD,size,ierr) |
257 |
|
|
CALL MPI_Comm_rank(MPI_COMM_WORLD,rank,ierr) |
258 |
|
|
npd=size-(NCPUs_PRNT+NCPUs_CHLD(1)) |
259 |
|
|
if(rank.lt.npd) color=0 |
260 |
|
|
CALL MPI_COMM_SPLIT (MPI_COMM_WORLD, color,0, |
261 |
|
|
& NEST_COMM,ierr) |
262 |
|
|
|
263 |
|
|
CALL MPI_Comm_size(NEST_COMM,isize,ierr) |
264 |
|
|
CALL MPI_Comm_rank(NEST_COMM,irank,ierr) |
265 |
|
|
|
266 |
|
|
C-------------------------------------------------------- |
267 |
|
|
C- change local dir to rank_N and open log file |
268 |
|
|
CALL SETDIR( rank ) |
269 |
|
|
WRITE(fNam,'(A,I4.4)') 'STDlog.', rank |
270 |
|
|
OPEN( iUnit, FILE=fNam, STATUS='unknown') |
271 |
|
|
mLoop = 0 |
272 |
|
|
C-------------------------------------------------------- |
273 |
|
|
C COMPUTE MASTER VALUES |
274 |
|
|
C-------------------------------------------------------- |
275 |
|
|
MSTR_DRV(1) = 0 |
276 |
|
|
|
277 |
|
|
MSTR_PRNT(1) = npd |
278 |
|
|
MSTR_CHLD(1) = NCPUs_PRNT + npd |
279 |
|
|
|
280 |
|
|
DO Count_Lev = 2, NST_LEV_TOT |
281 |
|
|
MSTR_DRV(Count_Lev) = MSTR_CHLD(Count_Lev-1) + |
282 |
|
|
& NCPUs_CHLD(Count_Lev - 1) |
283 |
|
|
|
284 |
|
|
MSTR_CHLD(Count_Lev) = MSTR_DRV(Count_Lev) + 1 |
285 |
|
|
MSTR_PRNT(Count_Lev) = MSTR_CHLD(Count_Lev-1) |
286 |
|
|
ENDDO |
287 |
|
|
|
288 |
|
|
vstart = 1+rank*(nSxP/MSTR_PRNT(1)) |
289 |
|
|
vstop = (1+rank)*(nSxP/MSTR_PRNT(1)) |
290 |
|
|
VCONT = (nSxP/npd)*(nSyP)*rank-1 |
291 |
|
|
VCONTP(rank) = VCONT |
292 |
|
|
|
293 |
|
|
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
294 |
|
|
C-- Print out nesting parameter: |
295 |
|
|
|
296 |
|
|
WRITE(iUnit,'(A)') '// ===================================' |
297 |
|
|
WRITE(iUnit,'(A)') '// NEST_DRIVER parameters :' |
298 |
|
|
WRITE(iUnit,'(A)') '// ===================================' |
299 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: rank =', rank, ' ; color =',color |
300 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: size =', size, ' ; npd =', npd |
301 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: irank =', irank, ' ; isize =', isize |
302 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: vstart =', vstart |
303 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: vstop =', vstop |
304 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: VCONTP =', VCONTP(rank) |
305 |
|
|
|
306 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: NST_LEV_TOT =', NST_LEV_TOT |
307 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: NST_LEV =', NST_LEV |
308 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: NCPUs_PRNT =', NCPUs_PRNT |
309 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: NCPUs_CHLD =', NCPUs_CHLD |
310 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: MSTR_DRV =', MSTR_DRV |
311 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: MSTR_PRNT =', MSTR_PRNT |
312 |
|
|
WRITE(iUnit,*) 'NEST_DRIVER: MSTR_CHLD =', MSTR_CHLD |
313 |
|
|
|
314 |
|
|
C-------------------------------------------------------- |
315 |
|
|
C COMPUTE DOMAIN DECOMPOSITION |
316 |
|
|
C-------------------------------------------------------- |
317 |
|
|
c if(rank.eq.0) then |
318 |
|
|
|
319 |
|
|
IM_C = int(NxC/nSxC) |
320 |
|
|
JM_C = int(NyC/nSyC) |
321 |
|
|
|
322 |
|
|
IM_P = int(NxP/nSxP) |
323 |
|
|
JM_P = int(NyP/nSyP) |
324 |
|
|
|
325 |
|
|
ICONT = 0 |
326 |
|
|
DO I = 1,nSxP |
327 |
|
|
DO J = 1,nSyP |
328 |
|
|
ICONT = ICONT + 1 |
329 |
|
|
IndI_P(ICONT) = IM_P*(I-1) |
330 |
|
|
IndJ_P(ICONT) = JM_P*(J-1) |
331 |
|
|
END DO |
332 |
|
|
END DO |
333 |
|
|
|
334 |
|
|
ICONT = 0 |
335 |
|
|
DO I = 1,nSxC |
336 |
|
|
DO J = 1,nSyC |
337 |
|
|
ICONT = ICONT + 1 |
338 |
|
|
IndI_C(ICONT) = IM_C*(I-1) |
339 |
|
|
IndJ_C(ICONT) = JM_C*(J-1) |
340 |
|
|
END DO |
341 |
|
|
END DO |
342 |
|
|
|
343 |
|
|
index = 6*JM_P*NrP*5 |
344 |
|
|
index_var3D = IM_P*JM_P*NrP |
345 |
|
|
index_var2D = IM_P*JM_P |
346 |
|
|
|
347 |
|
|
indexF = (JM_C+OLY+OLY)*NrC*2*5 |
348 |
|
|
index3F = IM_C*JM_C*NrC*n3dC |
349 |
|
|
index2F = IM_C*JM_C*4 |
350 |
|
|
|
351 |
|
|
ALLOCATE( globalPA(6,JM_P,NrP,5) ) |
352 |
|
|
ALLOCATE( localP3D_a(IM_P,JM_P,NrP) ) |
353 |
|
|
ALLOCATE( localP2D_a(IM_P,JM_P) ) |
354 |
|
|
|
355 |
|
|
ALLOCATE( VAR_C1((JM_C+OLY+OLY),NrC,2,5) ) |
356 |
|
|
ALLOCATE( globalC3D_a(IM_C,JM_C,NrC,n3dC) ) |
357 |
|
|
ALLOCATE( globalC2D_a(IM_C,JM_C,4)) |
358 |
|
|
|
359 |
|
|
IF ( rank.EQ.0 ) THEN |
360 |
|
|
C-------------------------------------------------------- |
361 |
|
|
C WARNING |
362 |
|
|
C-------------------------------------------------------- |
363 |
|
|
write(iUnit,*) '*************************************' |
364 |
|
|
write(iUnit,*) ' have you update geometric files?' |
365 |
|
|
write(iUnit,*) ' in ./CHILD e ./PARENT' |
366 |
|
|
write(iUnit,*) '*************************************' |
367 |
|
|
C-------------------------------------------------------- |
368 |
|
|
C PARENT MODEL |
369 |
|
|
C-------------------------------------------------------- |
370 |
|
|
write(iUnit,*) ' [1] Read PARENT model geometry' |
371 |
|
|
C---------------------------------------------------- |
372 |
|
|
C XC & YC |
373 |
|
|
C---------------------------------------------------- |
374 |
|
|
MSIZE = NxP*NyP*WORDLENGTH |
375 |
|
|
|
376 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
377 |
|
|
& file=trim(dirNEST)//'/PARENT/XC.data', |
378 |
|
|
& form='unformatted') |
379 |
|
|
|
380 |
|
|
read (1,REC=1) Xo_P(:,:) |
381 |
|
|
close(1) |
382 |
|
|
|
383 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
384 |
|
|
& file=trim(dirNEST)//'/PARENT/YC.data', |
385 |
|
|
& form='unformatted') |
386 |
|
|
|
387 |
|
|
read (1,REC=1) Yo_P(:,:) |
388 |
|
|
close(1) |
389 |
|
|
C---------------------------------------------------- |
390 |
|
|
C XG & YG |
391 |
|
|
C---------------------------------------------------- |
392 |
|
|
MSIZE = NxP*NyP*WORDLENGTH |
393 |
|
|
|
394 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
395 |
|
|
& file=trim(dirNEST)//'/PARENT/XG.data', |
396 |
|
|
& form='unformatted') |
397 |
|
|
|
398 |
|
|
read (1,REC=1) Xg_P(:,:) |
399 |
|
|
close(1) |
400 |
|
|
|
401 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
402 |
|
|
& file=trim(dirNEST)//'/PARENT/YG.data', |
403 |
|
|
& form='unformatted') |
404 |
|
|
|
405 |
|
|
read (1,REC=1) Yg_P(:,:) |
406 |
|
|
close(1) |
407 |
|
|
C---------------------------------------------------- |
408 |
|
|
C Yu |
409 |
|
|
C---------------------------------------------------- |
410 |
|
|
DO J = 1,NyP |
411 |
|
|
DO I = 1,NxP |
412 |
|
|
c Xu_P(I,J) = Xg_P(I,J) |
413 |
|
|
Yu_P(I,J) = Yo_P(I,J) |
414 |
|
|
ENDDO |
415 |
|
|
ENDDO |
416 |
|
|
C---------------------------------------------------- |
417 |
|
|
C Xv & Yv |
418 |
|
|
C---------------------------------------------------- |
419 |
|
|
DO J = 1,NyP |
420 |
|
|
DO I = 1,NxP |
421 |
|
|
Xv_P(I,J) = Xo_P(I,J) |
422 |
|
|
Yv_P(I,J) = Yg_P(I,J) |
423 |
|
|
ENDDO |
424 |
|
|
ENDDO |
425 |
|
|
C---------------------------------------------------- |
426 |
|
|
C hFacC |
427 |
|
|
C---------------------------------------------------- |
428 |
|
|
MSIZE = NxP*NyP*NrP*WORDLENGTH |
429 |
|
|
|
430 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
431 |
|
|
& file=trim(dirNEST)//'/PARENT/hFacC.data', |
432 |
|
|
& form='unformatted') |
433 |
|
|
|
434 |
|
|
read (1,REC=1) hFacC_P(:,:,:) |
435 |
|
|
close(1) |
436 |
|
|
C---------------------------------------------------- |
437 |
|
|
C hFacW |
438 |
|
|
C---------------------------------------------------- |
439 |
|
|
MSIZE = NxP*NyP*NrP*WORDLENGTH |
440 |
|
|
|
441 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
442 |
|
|
& file=trim(dirNEST)//'/PARENT/hFacW.data', |
443 |
|
|
& form='unformatted') |
444 |
|
|
|
445 |
|
|
read (1,REC=1) hFacW_P(:,:,:) |
446 |
|
|
close(1) |
447 |
|
|
C---------------------------------------------------- |
448 |
|
|
C hFacS |
449 |
|
|
C---------------------------------------------------- |
450 |
|
|
MSIZE = NxP*NyP*NrP*WORDLENGTH |
451 |
|
|
|
452 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
453 |
|
|
& file=trim(dirNEST)//'/PARENT/hFacS.data', |
454 |
|
|
& form='unformatted') |
455 |
|
|
|
456 |
|
|
read (1,REC=1) hFacS_P(:,:,:) |
457 |
|
|
close(1) |
458 |
|
|
C---------------------------------------------------- |
459 |
|
|
C RAC |
460 |
|
|
C---------------------------------------------------- |
461 |
|
|
MSIZE = NxP*NyP*WORDLENGTH |
462 |
|
|
|
463 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
464 |
|
|
& file=trim(dirNEST)//'/PARENT/RAC.data', |
465 |
|
|
& form='unformatted') |
466 |
|
|
|
467 |
|
|
read (1,REC=1) RAC_P(:,:) |
468 |
|
|
close(1) |
469 |
|
|
C---------------------------------------------------- |
470 |
|
|
C RAW |
471 |
|
|
C---------------------------------------------------- |
472 |
|
|
MSIZE = NxP*NyP*WORDLENGTH |
473 |
|
|
|
474 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
475 |
|
|
& file=trim(dirNEST)//'/PARENT/RAW.data', |
476 |
|
|
& form='unformatted') |
477 |
|
|
|
478 |
|
|
read (1,REC=1) RAW_P(:,:) |
479 |
|
|
close(1) |
480 |
|
|
C---------------------------------------------------- |
481 |
|
|
C RAS |
482 |
|
|
C---------------------------------------------------- |
483 |
|
|
MSIZE = NxP*NyP*WORDLENGTH |
484 |
|
|
|
485 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
486 |
|
|
& file=trim(dirNEST)//'/PARENT/RAS.data', |
487 |
|
|
& form='unformatted') |
488 |
|
|
|
489 |
|
|
read (1,REC=1) RAS_P(:,:) |
490 |
|
|
close(1) |
491 |
|
|
C---------------------------------------------------- |
492 |
|
|
C MASK x PARENT |
493 |
|
|
C---------------------------------------------------- |
494 |
|
|
DO K = 1,NrP |
495 |
|
|
DO J = 1,NyP |
496 |
|
|
DO I = 1,NxP |
497 |
|
|
DEEP_P(i,j,k) = 0. |
498 |
|
|
IF (hFacC_P(i,j,k).ne.0) then |
499 |
|
|
DEEP_P(I,J,K) = 1. |
500 |
|
|
ENDIF |
501 |
|
|
ENDDO |
502 |
|
|
ENDDO |
503 |
|
|
ENDDO |
504 |
|
|
C---------------------------------------------------- |
505 |
|
|
C 1/Volume (C) |
506 |
|
|
C---------------------------------------------------- |
507 |
|
|
DO K = 1,NrP |
508 |
|
|
DO J = 1,NyP |
509 |
|
|
DO I = 1,NxP |
510 |
|
|
INV_VOL_C_P(I,J,K) = 1. |
511 |
|
|
IF ((RAC_P(I,J)*hFacC_P(I,J,K)).ne.0.) THEN |
512 |
|
|
INV_VOL_C_P(I,J,K) = 1./(RAC_P(I,J)*hFacC_P(I,J,K)) |
513 |
|
|
ENDIF |
514 |
|
|
ENDDO |
515 |
|
|
ENDDO |
516 |
|
|
ENDDO |
517 |
|
|
C---------------------------------------------------- |
518 |
|
|
C 1/Volume (W) |
519 |
|
|
C---------------------------------------------------- |
520 |
|
|
DO K = 1,NrP |
521 |
|
|
DO J = 1,NyP |
522 |
|
|
DO I = 1,NxP |
523 |
|
|
INV_VOL_W_P(I,J,K) = 1. |
524 |
|
|
IF ((RAW_P(I,J)*hFacW_P(I,J,K)).ne.0.) THEN |
525 |
|
|
INV_VOL_W_P(I,J,K) = 1./(RAW_P(I,J)*hFacW_P(I,J,K)) |
526 |
|
|
ENDIF |
527 |
|
|
ENDDO |
528 |
|
|
ENDDO |
529 |
|
|
ENDDO |
530 |
|
|
C---------------------------------------------------- |
531 |
|
|
C 1/Volume (S) |
532 |
|
|
C---------------------------------------------------- |
533 |
|
|
DO K = 1,NrP |
534 |
|
|
DO J = 1,NyP |
535 |
|
|
DO I = 1,NxP |
536 |
|
|
INV_VOL_S_P(I,J,K) = 1. |
537 |
|
|
IF ((RAS_P(I,J)*hFacS_P(I,J,K)).ne.0.) THEN |
538 |
|
|
INV_VOL_S_P(I,J,K) = 1./(RAS_P(I,J)*hFacS_P(I,J,K)) |
539 |
|
|
ENDIF |
540 |
|
|
ENDDO |
541 |
|
|
ENDDO |
542 |
|
|
ENDDO |
543 |
|
|
C-------------------------------------------------------- |
544 |
|
|
C CHILD MODEL |
545 |
|
|
C-------------------------------------------------------- |
546 |
|
|
write(iUnit,*) ' [2] Read CHILD model geometry' |
547 |
|
|
C---------------------------------------------------- |
548 |
|
|
C XC & YC |
549 |
|
|
C---------------------------------------------------- |
550 |
|
|
MSIZE = NxC*NyC*WORDLENGTH |
551 |
|
|
|
552 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
553 |
|
|
& file=trim(dirNEST)//'/CHILD/XC.data', |
554 |
|
|
& form='unformatted') |
555 |
|
|
|
556 |
|
|
read (1,REC=1) Xo_C(:,:) |
557 |
|
|
close(1) |
558 |
|
|
|
559 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
560 |
|
|
& file=trim(dirNEST)//'/CHILD/YC.data', |
561 |
|
|
& form='unformatted') |
562 |
|
|
|
563 |
|
|
read (1,REC=1) Yo_C(:,:) |
564 |
|
|
close(1) |
565 |
|
|
C---------------------------------------------------- |
566 |
|
|
C XG & YG |
567 |
|
|
C---------------------------------------------------- |
568 |
|
|
MSIZE = NxC*NyC*WORDLENGTH |
569 |
|
|
|
570 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
571 |
|
|
& file=trim(dirNEST)//'/CHILD/XG.data', |
572 |
|
|
& form='unformatted') |
573 |
|
|
|
574 |
|
|
read (1,REC=1) Xg_C(:,:) |
575 |
|
|
close(1) |
576 |
|
|
|
577 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
578 |
|
|
& file=trim(dirNEST)//'/CHILD/YG.data', |
579 |
|
|
& form='unformatted') |
580 |
|
|
|
581 |
|
|
read (1,REC=1) Yg_C(:,:) |
582 |
|
|
close(1) |
583 |
|
|
C---------------------------------------------------- |
584 |
|
|
C Xu & Yu |
585 |
|
|
C---------------------------------------------------- |
586 |
|
|
DO J = 1,NyC |
587 |
|
|
DO I = 1,NxC |
588 |
|
|
Xu_C(I,J) = Xg_C(I,J) |
589 |
|
|
Yu_C(I,J) = Yo_C(I,J) |
590 |
|
|
ENDDO |
591 |
|
|
ENDDO |
592 |
|
|
C---------------------------------------------------- |
593 |
|
|
C Xv & Yv |
594 |
|
|
C---------------------------------------------------- |
595 |
|
|
DO J = 1,NyC |
596 |
|
|
DO I = 1,NxC |
597 |
|
|
Xv_C(I,J) = Xo_C(I,J) |
598 |
|
|
Yv_C(I,J) = Yg_C(I,J) |
599 |
|
|
ENDDO |
600 |
|
|
ENDDO |
601 |
|
|
C---------------------------------------------------- |
602 |
|
|
C hFacC |
603 |
|
|
C---------------------------------------------------- |
604 |
|
|
MSIZE = NxC*NyC*NrC*WORDLENGTH |
605 |
|
|
|
606 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
607 |
|
|
& file=trim(dirNEST)//'/CHILD/hFacC.data', |
608 |
|
|
& form='unformatted') |
609 |
|
|
|
610 |
|
|
read (1,REC=1) hFacC_C(:,:,:) |
611 |
|
|
close(1) |
612 |
|
|
C---------------------------------------------------- |
613 |
|
|
C hFacW |
614 |
|
|
C---------------------------------------------------- |
615 |
|
|
MSIZE = NxC*NyC*NrC*WORDLENGTH |
616 |
|
|
|
617 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
618 |
|
|
& file=trim(dirNEST)//'/CHILD/hFacW.data', |
619 |
|
|
& form='unformatted') |
620 |
|
|
|
621 |
|
|
read (1,REC=1) hFacW_C(:,:,:) |
622 |
|
|
close(1) |
623 |
|
|
C---------------------------------------------------- |
624 |
|
|
C hFacC |
625 |
|
|
C---------------------------------------------------- |
626 |
|
|
MSIZE = NxC*NyC*NrC*WORDLENGTH |
627 |
|
|
|
628 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
629 |
|
|
& file=trim(dirNEST)//'/CHILD/hFacS.data', |
630 |
|
|
& form='unformatted') |
631 |
|
|
|
632 |
|
|
read (1,REC=1) hFacS_C(:,:,:) |
633 |
|
|
close(1) |
634 |
|
|
C---------------------------------------------------- |
635 |
|
|
C RAC |
636 |
|
|
C---------------------------------------------------- |
637 |
|
|
MSIZE = NxC*NyC*WORDLENGTH |
638 |
|
|
|
639 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
640 |
|
|
& file=trim(dirNEST)//'/CHILD/RAC.data', |
641 |
|
|
& form='unformatted') |
642 |
|
|
|
643 |
|
|
read (1,REC=1) RAC_C(:,:) |
644 |
|
|
close(1) |
645 |
|
|
C---------------------------------------------------- |
646 |
|
|
C RAW |
647 |
|
|
C---------------------------------------------------- |
648 |
|
|
MSIZE = NxC*NyC*WORDLENGTH |
649 |
|
|
|
650 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
651 |
|
|
& file=trim(dirNEST)//'/CHILD/RAW.data', |
652 |
|
|
& form='unformatted') |
653 |
|
|
|
654 |
|
|
read (1,REC=1) RAW_C(:,:) |
655 |
|
|
close(1) |
656 |
|
|
C---------------------------------------------------- |
657 |
|
|
C RAS |
658 |
|
|
C---------------------------------------------------- |
659 |
|
|
MSIZE = NxC*NyC*WORDLENGTH |
660 |
|
|
|
661 |
|
|
open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD', |
662 |
|
|
& file=trim(dirNEST)//'/CHILD/RAS.data', |
663 |
|
|
& form='unformatted') |
664 |
|
|
|
665 |
|
|
read (1,REC=1) RAS_C(:,:) |
666 |
|
|
close(1) |
667 |
|
|
C---------------------------------------------------- |
668 |
|
|
C MASK x CHILD |
669 |
|
|
C---------------------------------------------------- |
670 |
|
|
DO K = 1,NrC |
671 |
|
|
DO J = 1,NyC |
672 |
|
|
DO I = 1,NxC |
673 |
|
|
DEEP_C(i,j,k) = 0. |
674 |
|
|
IF (hFacC_C(i,j,k).ne.0) then |
675 |
|
|
DEEP_C(I,J,K) = 1. |
676 |
|
|
ENDIF |
677 |
|
|
ENDDO |
678 |
|
|
ENDDO |
679 |
|
|
ENDDO |
680 |
|
|
|
681 |
|
|
C---------------------------------------------------- |
682 |
|
|
C __/________ ___________ |
683 |
|
|
C | | | | || |
684 |
|
|
C > o > | | | |
685 |
|
|
C |__/__|_____| |
686 |
|
|
C | | | |
687 |
|
|
C > o > | |
688 |
|
|
C |_____|_____|_____|_____| |
689 |
|
|
C |
690 |
|
|
C---------------------------------------------------- |
691 |
|
|
write(iUnit,*) ' [3] Compute J index P-->C' |
692 |
|
|
C-------------------------------------------------------- |
693 |
|
|
C Compute J indicies for mapping P->C (I) |
694 |
|
|
C-------------------------------------------------------- |
695 |
|
|
I = 1 |
696 |
|
|
II = WesternB |
697 |
|
|
|
698 |
|
|
DO J = 1,NyC |
699 |
|
|
P2C_U(J) = 0. |
700 |
|
|
DO JJ = 1,NyP-1 |
701 |
|
|
YF = Yg_C(I,J) |
702 |
|
|
YP1 = Yg_P(II,JJ) |
703 |
|
|
YP2 = Yg_P(II,JJ+1) |
704 |
|
|
IF (YF.ge.YP1.and.YF.lt.YP2) THEN |
705 |
|
|
P2C_U(J) = JJ |
706 |
|
|
ENDIF |
707 |
|
|
ENDDO |
708 |
|
|
ENDDO |
709 |
|
|
C-------------------------------------------------------- |
710 |
|
|
C Compute J indicies for mapping P->C (II) |
711 |
|
|
C-------------------------------------------------------- |
712 |
|
|
I = 1 |
713 |
|
|
II = WesternB |
714 |
|
|
|
715 |
|
|
DO J = 1,NyC |
716 |
|
|
P2C_linU(J) = 0. |
717 |
|
|
DO JJ = 1,NyP-1 |
718 |
|
|
YF = Yu_C(I,J) |
719 |
|
|
YP1 = Yu_P(II,JJ) |
720 |
|
|
YP2 = Yu_P(II,JJ+1) |
721 |
|
|
IF (YF.ge.YP1.and.YF.lt.YP2) THEN |
722 |
|
|
P2C_linU(J) = JJ |
723 |
|
|
ENDIF |
724 |
|
|
ENDDO |
725 |
|
|
ENDDO |
726 |
|
|
C-------------------------------------------------------- |
727 |
|
|
C Compute J indicies for mapping P->C (III) |
728 |
|
|
C-------------------------------------------------------- |
729 |
|
|
I = 1 |
730 |
|
|
II = WesternB |
731 |
|
|
|
732 |
|
|
DO J = 1,NyC |
733 |
|
|
DO JJ = 1,NyP-1 |
734 |
|
|
YF = Yu_C(I,J) |
735 |
|
|
YP1 = Yu_P(II,JJ) |
736 |
|
|
IF (YF.eq.YP1) THEN |
737 |
|
|
WO3_linU(J) = 0 |
738 |
|
|
if (J+1.le.NyC) WO3_linU(J+1) = 1 |
739 |
|
|
if (J+2.le.NyC) WO3_linU(J+2) = 2 |
740 |
|
|
ENDIF |
741 |
|
|
ENDDO |
742 |
|
|
ENDDO |
743 |
|
|
C--------------------Lower bound |
744 |
|
|
DO J = 1,NyC |
745 |
|
|
DO JJ = 1,NyP-1 |
746 |
|
|
YF = Yu_C(I,J) |
747 |
|
|
YP1 = Yu_P(II,JJ) |
748 |
|
|
IF (YF.eq.YP1) THEN |
749 |
|
|
WO3_linU(J) = 0 |
750 |
|
|
if (J-1.gt.0) WO3_linU(J-1) = 2 |
751 |
|
|
if (J-2.gt.0) WO3_linU(J-2) = 1 |
752 |
|
|
GOTO 2345 |
753 |
|
|
ENDIF |
754 |
|
|
ENDDO |
755 |
|
|
ENDDO |
756 |
|
|
2345 CONTINUE |
757 |
|
|
C--------------------Upper bound |
758 |
|
|
DO J = NyC,1,-1 |
759 |
|
|
DO JJ = 1,NyP-1 |
760 |
|
|
YF = Yu_C(I,J) |
761 |
|
|
YP1 = Yu_P(II,JJ) |
762 |
|
|
IF (YF.eq.YP1) THEN |
763 |
|
|
WO3_linU(J) = 0 |
764 |
|
|
if (J+1.le.NyC) WO3_linU(J+1) = 1 |
765 |
|
|
if (J+2.le.NyC) WO3_linU(J+2) = 2 |
766 |
|
|
GOTO 2346 |
767 |
|
|
ENDIF |
768 |
|
|
ENDDO |
769 |
|
|
ENDDO |
770 |
|
|
2346 CONTINUE |
771 |
|
|
C-------------------------------------------------------- |
772 |
|
|
C Compute J indicies for mapping P->C (IV) |
773 |
|
|
C-------------------------------------------------------- |
774 |
|
|
I = 1 |
775 |
|
|
II = WesternB |
776 |
|
|
|
777 |
|
|
DO J = 1,NyC |
778 |
|
|
P2C_linV(J) = 0. |
779 |
|
|
DO JJ = 1,NyP-1 |
780 |
|
|
YF = Yv_C(I,J) |
781 |
|
|
YP1 = Yv_P(II,JJ) |
782 |
|
|
YP2 = Yv_P(II,JJ+1) |
783 |
|
|
IF (YF.ge.YP1.and.YF.lt.YP2) THEN |
784 |
|
|
P2C_linV(J) = JJ |
785 |
|
|
ENDIF |
786 |
|
|
ENDDO |
787 |
|
|
ENDDO |
788 |
|
|
C-------------------------------------------------------- |
789 |
|
|
C Compute J indicies for mapping P->C (V) |
790 |
|
|
C-------------------------------------------------------- |
791 |
|
|
I = 1 |
792 |
|
|
II = WesternB |
793 |
|
|
|
794 |
|
|
DO J = 1,NyC |
795 |
|
|
DO JJ = 1,NyP-1 |
796 |
|
|
YF = Yv_C(I,J) |
797 |
|
|
YP1 = Yv_P(II,JJ) |
798 |
|
|
IF (YF.eq.YP1) THEN |
799 |
|
|
WO3_linV(J) = 0 |
800 |
|
|
if (J+1.le.NyC) WO3_linV(J+1) = 1 |
801 |
|
|
if (J+2.le.NyC) WO3_linV(J+2) = 2 |
802 |
|
|
ENDIF |
803 |
|
|
ENDDO |
804 |
|
|
ENDDO |
805 |
|
|
C--------------------Lower bound |
806 |
|
|
DO J = 1,NyC |
807 |
|
|
DO JJ = 1,NyP-1 |
808 |
|
|
YF = Yv_C(I,J) |
809 |
|
|
YP1 = Yv_P(II,JJ) |
810 |
|
|
IF (YF.eq.YP1) THEN |
811 |
|
|
WO3_linV(J) = 0 |
812 |
|
|
if (J-1.gt.0) WO3_linV(J-1) = 2 |
813 |
|
|
if (J-2.gt.0) WO3_linV(J-2) = 1 |
814 |
|
|
GOTO 23451 |
815 |
|
|
ENDIF |
816 |
|
|
ENDDO |
817 |
|
|
ENDDO |
818 |
|
|
23451 CONTINUE |
819 |
|
|
C--------------------Upper bound |
820 |
|
|
DO J = NyC,1,-1 |
821 |
|
|
DO JJ = 1,NyP-1 |
822 |
|
|
YF = Yv_C(I,J) |
823 |
|
|
YP1 = Yv_P(II,JJ) |
824 |
|
|
IF (YF.eq.YP1) THEN |
825 |
|
|
WO3_linV(J) = 0 |
826 |
|
|
if (J+1.le.NyC) WO3_linV(J+1) = 1 |
827 |
|
|
if (J+2.le.NyC) WO3_linV(J+2) = 2 |
828 |
|
|
GOTO 23461 |
829 |
|
|
ENDIF |
830 |
|
|
ENDDO |
831 |
|
|
ENDDO |
832 |
|
|
23461 CONTINUE |
833 |
|
|
C-------------------------------------------------------- |
834 |
|
|
C Compute J indicies for mapping P->C (V) |
835 |
|
|
C-------------------------------------------------------- |
836 |
|
|
write(iUnit,*) ' [5] Compute J index P-->C for (o)' |
837 |
|
|
I = 1 |
838 |
|
|
II = WesternB |
839 |
|
|
|
840 |
|
|
DO J = 1,NyC |
841 |
|
|
P2C_o(J) = 0. |
842 |
|
|
DO JJ = 1,NyP-1 |
843 |
|
|
YF = Yo_C(I,J) |
844 |
|
|
YP1 = Yg_P(II,JJ) |
845 |
|
|
YP2 = Yg_P(II,JJ+1) |
846 |
|
|
IF (YF.gt.YP1.and.YF.lt.YP2) THEN |
847 |
|
|
P2C_o(J) = JJ |
848 |
|
|
ENDIF |
849 |
|
|
ENDDO |
850 |
|
|
ENDDO |
851 |
|
|
C-------------------------------------------------------- |
852 |
|
|
C Compute J indicies for mapping P->C (VI) |
853 |
|
|
C-------------------------------------------------------- |
854 |
|
|
write(iUnit,*) ' [6] Compute J index P-->C for (v bilinear)' |
855 |
|
|
I = 1 |
856 |
|
|
II = WesternB |
857 |
|
|
|
858 |
|
|
DO J = 1,NyC |
859 |
|
|
DO JJ = 2,NyP-1 |
860 |
|
|
YF = Yv_C(I,J) |
861 |
|
|
YP1 = Yv_P(II,JJ) |
862 |
|
|
YP2 = Yv_P(II,JJ+1) |
863 |
|
|
YP3 = Yv_P(II,JJ-1) |
864 |
|
|
|
865 |
|
|
IF (YF.ge.YP1.and.YF.lt.YP2) THEN |
866 |
|
|
P2C1_V(J) = JJ |
867 |
|
|
P2C2_V(J) = JJ+1 |
868 |
|
|
ENDIF |
869 |
|
|
ENDDO |
870 |
|
|
ENDDO |
871 |
|
|
C-------------------------------------------------------- |
872 |
|
|
C Look for the 9 CHILD indicies in PARENT grid cell |
873 |
|
|
C-------------------------------------------------------- |
874 |
|
|
write(iUnit,*) ' [8] Compute I J index C-->P for (o)' |
875 |
|
|
|
876 |
|
|
DO J = 1,NyP |
877 |
|
|
DO I = 1,NxP |
878 |
|
|
I_C2P(:,I,J) = 0 |
879 |
|
|
J_C2P(:,I,J) = 0 |
880 |
|
|
|
881 |
|
|
DO JJ = 1,NyC |
882 |
|
|
DO II = 1,NxC |
883 |
|
|
IF (Xo_C(II,JJ).eq.Xo_P(I,J).and. |
884 |
|
|
& Yo_C(II,JJ).eq.Yo_P(I,J)) then |
885 |
|
|
|
886 |
|
|
KK = 0 |
887 |
|
|
DO JJJ = JJ-1,JJ+1 |
888 |
|
|
DO III = II-1,II+1 |
889 |
|
|
KK = kk +1 |
890 |
|
|
if (III.lt.1.or.III.gt.NxC) cycle |
891 |
|
|
if (JJJ.lt.1.or.JJJ.gt.NyC) cycle |
892 |
|
|
I_C2P(KK,I,J) = III |
893 |
|
|
J_C2P(KK,I,J) = JJJ |
894 |
|
|
ENDDO |
895 |
|
|
ENDDO |
896 |
|
|
ENDIF |
897 |
|
|
|
898 |
|
|
ENDDO |
899 |
|
|
ENDDO |
900 |
|
|
|
901 |
|
|
ENDDO |
902 |
|
|
ENDDO |
903 |
|
|
C-- end if rank=0 |
904 |
|
|
ENDIF |
905 |
|
|
C-------------------------------------------------------- |
906 |
|
|
C Broadcast all the above variables |
907 |
|
|
C-------------------------------------------------------- |
908 |
|
|
CALL MPI_BCAST(I_C2P,9*NxP*NyP,MPI_INTEGER, |
909 |
|
|
& 0,NEST_COMM,ierr) |
910 |
|
|
CALL MPI_BCAST(J_C2P,9*NxP*NyP,MPI_INTEGER, |
911 |
|
|
& 0,NEST_COMM,ierr) |
912 |
|
|
|
913 |
|
|
CALL MPI_BCAST(RAC_C,NxC*NyC,MPI_REAL, |
914 |
|
|
& 0,NEST_COMM,ierr) |
915 |
|
|
CALL MPI_BCAST(hFacC_C,NxC*NyC*NrC,MPI_REAL, |
916 |
|
|
& 0,NEST_COMM,ierr) |
917 |
|
|
CALL MPI_BCAST(INV_VOL_C_P,NxP*NyP*NrP,MPI_REAL, |
918 |
|
|
& 0,NEST_COMM,ierr) |
919 |
|
|
|
920 |
|
|
CALL MPI_BCAST(RAW_C,NxC*NyC,MPI_REAL, |
921 |
|
|
& 0,NEST_COMM,ierr) |
922 |
|
|
CALL MPI_BCAST(hFacW_C,NxC*NyC*NrC,MPI_REAL, |
923 |
|
|
& 0,NEST_COMM,ierr) |
924 |
|
|
CALL MPI_BCAST(INV_VOL_W_P,NxP*NyP*NrP,MPI_REAL, |
925 |
|
|
& 0,NEST_COMM,ierr) |
926 |
|
|
|
927 |
|
|
CALL MPI_BCAST(RAS_C,NxC*NyC,MPI_REAL, |
928 |
|
|
& 0,NEST_COMM,ierr) |
929 |
|
|
CALL MPI_BCAST(hFacS_C,NxC*NyC*NrC,MPI_REAL, |
930 |
|
|
& 0,NEST_COMM,ierr) |
931 |
|
|
CALL MPI_BCAST(INV_VOL_S_P,NxP*NyP*NrP,MPI_REAL, |
932 |
|
|
& 0,NEST_COMM,ierr) |
933 |
|
|
|
934 |
|
|
CALL MPI_BCAST(DEEP_C,NxC*NyC*NrC,MPI_REAL, |
935 |
|
|
& 0,NEST_COMM,ierr) |
936 |
|
|
CALL MPI_BCAST(RAC_P,NxP*NyP,MPI_REAL, |
937 |
|
|
& 0,NEST_COMM,ierr) |
938 |
|
|
|
939 |
|
|
CALL MPI_BCAST(IM_P,1,MPI_INTEGER, |
940 |
|
|
& 0,NEST_COMM,ierr) |
941 |
|
|
CALL MPI_BCAST(JM_P,1,MPI_INTEGER, |
942 |
|
|
& 0,NEST_COMM,ierr) |
943 |
|
|
CALL MPI_BCAST(index_var3D,1,MPI_INTEGER, |
944 |
|
|
& 0,NEST_COMM,ierr) |
945 |
|
|
CALL MPI_BCAST(index_var2D,1,MPI_INTEGER, |
946 |
|
|
& 0,NEST_COMM,ierr) |
947 |
|
|
|
948 |
|
|
CALL MPI_BCAST(DEEP_P,NxP*NyP*NrP,MPI_REAL, |
949 |
|
|
& 0,NEST_COMM,ierr) |
950 |
|
|
CALL MPI_BCAST(hFacS_P,NxP*NyP*NrP,MPI_REAL, |
951 |
|
|
& 0,NEST_COMM,ierr) |
952 |
|
|
CALL MPI_BCAST(hFacC_P,NxP*NyP*NrP,MPI_REAL, |
953 |
|
|
& 0,NEST_COMM,ierr) |
954 |
|
|
CALL MPI_BCAST(hFacW_P,NxP*NyP*NrP,MPI_REAL, |
955 |
|
|
& 0,NEST_COMM,ierr) |
956 |
|
|
|
957 |
|
|
C-------------------------------------------------------- |
958 |
|
|
if(rank.eq.0) then |
959 |
|
|
C-------------------------------------------------------- |
960 |
|
|
DO K = 1,NrP |
961 |
|
|
DO J = 1,NyP |
962 |
|
|
DO I = WesternB+1,EasternB-1 |
963 |
|
|
C- WesternB side |
964 |
|
|
|
965 |
|
|
DO II = 1,9 |
966 |
|
|
IF (I_C2P(II,I,J).eq.0) cycle |
967 |
|
|
IF (J_C2P(II,I,J).eq.0) cycle |
968 |
|
|
|
969 |
|
|
Indx = I_C2P(II,I,J) |
970 |
|
|
Jndx = J_C2P(II,I,J) |
971 |
|
|
ENDDO |
972 |
|
|
ENDDO |
973 |
|
|
ENDDO |
974 |
|
|
ENDDO |
975 |
|
|
C--------------------------------------------------------- |
976 |
|
|
ONOFF=0 |
977 |
|
|
endif |
978 |
|
|
C-------------------------------------------------------- |
979 |
|
|
C BEGIN MAIN LOOP |
980 |
|
|
C-------------------------------------------------------- |
981 |
|
|
DO |
982 |
|
|
mLoop = mLoop + 1 |
983 |
|
|
c DO mLoop=1,24 |
984 |
|
|
WRITE(iUnit,'(A,I6)') '== Main Loop , iter=', mLoop |
985 |
|
|
if(rank.eq.0) then |
986 |
|
|
C-------------------------------------------------------- |
987 |
|
|
C (1) READ FROM PARENT MODEL |
988 |
|
|
C-------------------------------------------------------- |
989 |
|
|
ICONT=1 |
990 |
|
|
DO WHILE(ICONT.le.nSxP*nSyP) |
991 |
|
|
from= MPI_ANY_SOURCE |
992 |
|
|
|
993 |
|
|
CALL MPI_RECV (globalPA, index, MPI_REAL8, |
994 |
|
|
& FROM, 3000, |
995 |
|
|
& MPI_COMM_World, status,ierr) |
996 |
|
|
|
997 |
|
|
ICONT=ICONT+1 |
998 |
|
|
|
999 |
|
|
whm=status(MPI_SOURCE)-MSTR_PRNT(NST_LEV)+1 |
1000 |
|
|
|
1001 |
|
|
CALL MPI_GET_COUNT(status,MPI_REAL8,st_count,ierr) |
1002 |
|
|
|
1003 |
|
|
DO II = 1,6 |
1004 |
|
|
IF (globalPA(II,1,1,1).ne.-999.) THEN |
1005 |
|
|
globalP1(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1006 |
|
|
& globalPA(II,1:JM_P,:,1) |
1007 |
|
|
globalP2(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1008 |
|
|
& globalPA(II,1:JM_P,:,2) |
1009 |
|
|
globalP3(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1010 |
|
|
& globalPA(II,1:JM_P,:,3) |
1011 |
|
|
globalP4(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1012 |
|
|
& globalPA(II,1:JM_P,:,4) |
1013 |
|
|
globalP5(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1014 |
|
|
& globalPA(II,1:JM_P,:,5) |
1015 |
|
|
ENDIF |
1016 |
|
|
|
1017 |
|
|
ENDDO |
1018 |
|
|
ENDDO |
1019 |
|
|
C-------------------------------------------------------- |
1020 |
|
|
C Start interpolation for CHILD |
1021 |
|
|
C-------------------------------------------------------- |
1022 |
|
|
CALL INTERPOLATION_P2C ( |
1023 |
|
|
& globalP1,globalP2,globalP3,globalP4,globalP5, |
1024 |
|
|
& NxP,NyP,NrP, |
1025 |
|
|
& NxC,NyC,NrC, |
1026 |
|
|
$ WesternB,EasternB, |
1027 |
|
|
$ P2C_U,P2C_V,P2C_o,P2C1_V,P2C2_V,P2C1_o,P2C2_o, |
1028 |
|
|
$ P2C_linU,WO3_linU,P2C_linV,WO3_linV, |
1029 |
|
|
$ Xv_C,Yv_C,Xv_P,Yv_P, |
1030 |
|
|
$ T_C1,S_C1,U_C1,V_C1,ETA_C1, |
1031 |
|
|
$ DEEP_C,DEEP_P |
1032 |
|
|
& ) |
1033 |
|
|
C============================================================== |
1034 |
|
|
C Open Files from PARENT MODEL |
1035 |
|
|
C============================================================== |
1036 |
|
|
ICONT=1 |
1037 |
|
|
|
1038 |
|
|
do while(ICONT.le.nSxP*nSyP) |
1039 |
|
|
from= MPI_ANY_SOURCE |
1040 |
|
|
|
1041 |
|
|
CALL MPI_RECV (globalPA, index, MPI_REAL8, |
1042 |
|
|
& FROM, 3000, |
1043 |
|
|
& MPI_COMM_World, status,ierr) |
1044 |
|
|
|
1045 |
|
|
ICONT=ICONT+1 |
1046 |
|
|
|
1047 |
|
|
whm=status(MPI_SOURCE)-MSTR_PRNT(NST_LEV)+1 |
1048 |
|
|
|
1049 |
|
|
CALL MPI_GET_COUNT(status,MPI_REAL8,st_count,ierr) |
1050 |
|
|
|
1051 |
|
|
DO II = 1,6 |
1052 |
|
|
IF (globalPA(II,1,1,1).ne.-999.) THEN |
1053 |
|
|
globalP1(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1054 |
|
|
& globalPA(II,1:JM_P,:,1) |
1055 |
|
|
globalP2(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1056 |
|
|
& globalPA(II,1:JM_P,:,2) |
1057 |
|
|
globalP3(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1058 |
|
|
& globalPA(II,1:JM_P,:,3) |
1059 |
|
|
globalP4(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1060 |
|
|
& globalPA(II,1:JM_P,:,4) |
1061 |
|
|
globalP5(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) = |
1062 |
|
|
& globalPA(II,1:JM_P,:,5) |
1063 |
|
|
ENDIF |
1064 |
|
|
ENDDO |
1065 |
|
|
|
1066 |
|
|
end do |
1067 |
|
|
C-------------------------------------------------------- |
1068 |
|
|
C Start inteprolation for CHILD |
1069 |
|
|
C-------------------------------------------------------- |
1070 |
|
|
CALL INTERPOLATION_P2C ( |
1071 |
|
|
& globalP1,globalP2,globalP3,globalP4,globalP5, |
1072 |
|
|
& NxP,NyP,NrP, |
1073 |
|
|
& NxC,NyC,NrC, |
1074 |
|
|
$ WesternB,EasternB, |
1075 |
|
|
$ P2C_U,P2C_V,P2C_o,P2C1_V,P2C2_V,P2C1_o,P2C2_o, |
1076 |
|
|
$ P2C_linU,WO3_linU,P2C_linV,WO3_linV, |
1077 |
|
|
$ Xv_C,Yv_C,Xv_P,Yv_P, |
1078 |
|
|
$ T_C2,S_C2,U_C2,V_C2,ETA_C2, |
1079 |
|
|
$ DEEP_C,DEEP_P |
1080 |
|
|
& ) |
1081 |
|
|
|
1082 |
|
|
C============================================================== |
1083 |
|
|
C Temporal Interpolation OBCs x CHILD MODEL |
1084 |
|
|
C============================================================== |
1085 |
|
|
C 0 1200 |
1086 |
|
|
C ---+--.--.--+---- Parent |
1087 |
|
|
C |
1088 |
|
|
C |--|--|-- |
1089 |
|
|
C 0 800 |
1090 |
|
|
C 400 |
1091 |
|
|
C------------------------------------------------------------ |
1092 |
|
|
DO I = 1,2 ! WesternB & EasternB |
1093 |
|
|
DIFF_T(:,:,I) = (T_C2(:,:,I) - T_C1(:,:,I))/3. |
1094 |
|
|
DIFF_S(:,:,I) = (S_C2(:,:,I) - S_C1(:,:,I))/3. |
1095 |
|
|
DIFF_U(:,:,I) = (U_C2(:,:,I) - U_C1(:,:,I))/3. |
1096 |
|
|
DIFF_V(:,:,I) = (V_C2(:,:,I) - V_C1(:,:,I))/3. |
1097 |
|
|
DIFF_eta(:,:,I) = (eta_C2(:,:,I) - eta_C1(:,:,I))/3. |
1098 |
|
|
ENDDO |
1099 |
|
|
C------------------------------------------------------------- |
1100 |
|
|
C Step 0 (Rec = 1 ==> WesternB) |
1101 |
|
|
C------- (Rec = 2 ==> EasternB) |
1102 |
|
|
|
1103 |
|
|
DO I = 1,2 !WesternB & EasternB |
1104 |
|
|
T_C1(:,:,I) = T_C2(:,:,I) !+ DIFF_T(:,:,I) |
1105 |
|
|
S_C1(:,:,I) = S_C2(:,:,I) !+ DIFF_S(:,:,I) |
1106 |
|
|
U_C1(:,:,I) = U_C2(:,:,I) !+ DIFF_U(:,:,I) |
1107 |
|
|
V_C1(:,:,I) = V_C2(:,:,I) !+ DIFF_V(:,:,I) |
1108 |
|
|
ETA_C1(:,:,I) = ETA_C2(:,:,I) !+ DIFF_ETA(:,:,I) |
1109 |
|
|
ENDDO |
1110 |
|
|
|
1111 |
|
|
if(ONOFF.eq.0) then |
1112 |
|
|
C--------------------------------------------------------------------- |
1113 |
|
|
ICONT = -1 |
1114 |
|
|
DO I = 1,nSxC |
1115 |
|
|
DO J = 1,nSyC |
1116 |
|
|
ICONT = ICONT + 1 |
1117 |
|
|
IndI = IM_C*(I-1) |
1118 |
|
|
IndJ = JM_C*(J-1) |
1119 |
|
|
|
1120 |
|
|
VAR_C1(:,:,:,:) = 0. |
1121 |
|
|
|
1122 |
|
|
J1 = 1+IndJ-OLY |
1123 |
|
|
J2 = JM_C+IndJ+OLY |
1124 |
|
|
|
1125 |
|
|
JJ1 = 1 |
1126 |
|
|
JJ2 = JM_C+OLY+OLY |
1127 |
|
|
|
1128 |
|
|
IF(1 +IndJ-OLY.LT.0) THEN |
1129 |
|
|
J1 = 1 |
1130 |
|
|
JJ1 = 4 |
1131 |
|
|
ENDIF |
1132 |
|
|
|
1133 |
|
|
IF(JM_C+IndJ+OLY.GT.NyC) THEN |
1134 |
|
|
J2 = NyC |
1135 |
|
|
JJ2 = JM_C |
1136 |
|
|
ENDIF |
1137 |
|
|
|
1138 |
|
|
VAR_C1(JJ1:JJ2,:,:,1) = U_C1(J1:J2,:,:) |
1139 |
|
|
VAR_C1(JJ1:JJ2,:,:,2) = V_C1(J1:J2,:,:) |
1140 |
|
|
VAR_C1(JJ1:JJ2,:,:,3) = T_C1(J1:J2,:,:) |
1141 |
|
|
VAR_C1(JJ1:JJ2,:,:,4) = S_C1(J1:J2,:,:) |
1142 |
|
|
VAR_C1(JJ1:JJ2,:,:,5) = ETA_C1(J1:J2,:,:) |
1143 |
|
|
|
1144 |
|
|
CALL MPI_SEND (VAR_C1, indexF, MPI_REAL8, |
1145 |
|
|
& MSTR_CHLD(NST_LEV)+ICONT, 3000, |
1146 |
|
|
& MPI_Comm_World,ierr) |
1147 |
|
|
|
1148 |
|
|
ENDDO |
1149 |
|
|
ENDDO |
1150 |
|
|
C---------------------------------------------------------------------- |
1151 |
|
|
c write(*,*) 'VIC: MANDO SEGNALE DI OK AL CHILD PER INIZIALIZZARE' |
1152 |
|
|
ONOFF=1 |
1153 |
|
|
ENDIF |
1154 |
|
|
c write(*,*) 'VIC: MANDO SEGNALE DI OK AL CHILD PER IL PASSO 1' |
1155 |
|
|
C----------------------------------------------------------------------- |
1156 |
|
|
ICONT = -1 |
1157 |
|
|
DO I = 1,nSxC |
1158 |
|
|
DO J = 1,nSyC |
1159 |
|
|
ICONT = ICONT + 1 |
1160 |
|
|
IndI = IM_C*(I-1) |
1161 |
|
|
IndJ = JM_C*(J-1) |
1162 |
|
|
|
1163 |
|
|
VAR_C1(:,:,:,:) = 0. |
1164 |
|
|
|
1165 |
|
|
J1 = 1+IndJ-OLY |
1166 |
|
|
J2 = JM_C+IndJ+OLY |
1167 |
|
|
|
1168 |
|
|
JJ1 = 1 |
1169 |
|
|
JJ2 = JM_C+OLY+OLY |
1170 |
|
|
|
1171 |
|
|
IF(1 +IndJ-OLY.LT.0) THEN |
1172 |
|
|
J1 = 1 |
1173 |
|
|
JJ1 = 4 |
1174 |
|
|
ENDIF |
1175 |
|
|
|
1176 |
|
|
IF(JM_C+IndJ+OLY.GT.NyC) THEN |
1177 |
|
|
J2 = NyC |
1178 |
|
|
JJ2 = JM_C |
1179 |
|
|
ENDIF |
1180 |
|
|
|
1181 |
|
|
VAR_C1(JJ1:JJ2,:,:,1) = U_C1(J1:J2,:,:) |
1182 |
|
|
VAR_C1(JJ1:JJ2,:,:,2) = V_C1(J1:J2,:,:) |
1183 |
|
|
VAR_C1(JJ1:JJ2,:,:,3) = T_C1(J1:J2,:,:) |
1184 |
|
|
VAR_C1(JJ1:JJ2,:,:,4) = S_C1(J1:J2,:,:) |
1185 |
|
|
VAR_C1(JJ1:JJ2,:,:,5) = ETA_C1(J1:J2,:,:) |
1186 |
|
|
|
1187 |
|
|
CALL MPI_SEND (VAR_C1, indexF, MPI_REAL8, |
1188 |
|
|
& MSTR_CHLD(NST_LEV)+ICONT, 3000, |
1189 |
|
|
& MPI_Comm_World,ierr) |
1190 |
|
|
|
1191 |
|
|
ENDDO |
1192 |
|
|
ENDDO |
1193 |
|
|
C--------------------------------------------------------------------- |
1194 |
|
|
goto 8888 |
1195 |
|
|
|
1196 |
|
|
C Step 1 (Rec = 3 ==> WesternB) |
1197 |
|
|
C------- (Rec = 4 ==> EasternB) |
1198 |
|
|
|
1199 |
|
|
DO I = 1,2 !WesternB & EasternB |
1200 |
|
|
T_C1(:,:,I) = T_C2(:,:,I) !+ DIFF_T(:,:,I) |
1201 |
|
|
S_C1(:,:,I) = S_C2(:,:,I) !+ DIFF_S(:,:,I) |
1202 |
|
|
U_C1(:,:,I) = U_C2(:,:,I) !+ DIFF_U(:,:,I) |
1203 |
|
|
V_C1(:,:,I) = V_C2(:,:,I) !+ DIFF_V(:,:,I) |
1204 |
|
|
ETA_C1(:,:,I) = ETA_C2(:,:,I) !+ DIFF_ETA(:,:,I) |
1205 |
|
|
ENDDO |
1206 |
|
|
C---------------------------------------------------------- |
1207 |
|
|
ICONT = -1 |
1208 |
|
|
DO I = 1,nSxC |
1209 |
|
|
DO J = 1,nSyC |
1210 |
|
|
ICONT = ICONT + 1 |
1211 |
|
|
IndI = IM_C*(I-1) |
1212 |
|
|
IndJ = JM_C*(J-1) |
1213 |
|
|
|
1214 |
|
|
VAR_C1(:,:,:,:) = 0. |
1215 |
|
|
|
1216 |
|
|
J1 = 1+IndJ-OLY |
1217 |
|
|
J2 = JM_C+IndJ+OLY |
1218 |
|
|
|
1219 |
|
|
JJ1 = 1 |
1220 |
|
|
JJ2 = JM_C+OLY+OLY |
1221 |
|
|
|
1222 |
|
|
IF(1 +IndJ-OLY.LT.0) THEN |
1223 |
|
|
J1 = 1 |
1224 |
|
|
JJ1 = 4 |
1225 |
|
|
ENDIF |
1226 |
|
|
|
1227 |
|
|
IF(JM_C+IndJ+OLY.GT.NyC) THEN |
1228 |
|
|
J2 = NyC |
1229 |
|
|
JJ2 = JM_C |
1230 |
|
|
ENDIF |
1231 |
|
|
|
1232 |
|
|
VAR_C1(JJ1:JJ2,:,:,1) = U_C1(J1:J2,:,:) |
1233 |
|
|
VAR_C1(JJ1:JJ2,:,:,2) = V_C1(J1:J2,:,:) |
1234 |
|
|
VAR_C1(JJ1:JJ2,:,:,3) = T_C1(J1:J2,:,:) |
1235 |
|
|
VAR_C1(JJ1:JJ2,:,:,4) = S_C1(J1:J2,:,:) |
1236 |
|
|
VAR_C1(JJ1:JJ2,:,:,5) = ETA_C1(J1:J2,:,:) |
1237 |
|
|
|
1238 |
|
|
CALL MPI_SEND (VAR_C1, indexF, MPI_REAL8, |
1239 |
|
|
& MSTR_CHLD(NST_LEV)+ICONT, 3000, |
1240 |
|
|
& MPI_Comm_World,ierr) |
1241 |
|
|
|
1242 |
|
|
ENDDO |
1243 |
|
|
ENDDO |
1244 |
|
|
C----------------------------------------------------------- |
1245 |
|
|
C Step 2 (Rec = 5 ==> WesternB) |
1246 |
|
|
C------- (Rec = 6 ==> EasternB) |
1247 |
|
|
|
1248 |
|
|
DO I = 1,2 !WesternB & EasternB |
1249 |
|
|
T_C1(:,:,I) = T_C2(:,:,I) !+ DIFF_T(:,:,I) |
1250 |
|
|
S_C1(:,:,I) = S_C2(:,:,I) !+ DIFF_S(:,:,I) |
1251 |
|
|
U_C1(:,:,I) = U_C2(:,:,I) !+ DIFF_U(:,:,I) |
1252 |
|
|
V_C1(:,:,I) = V_C2(:,:,I) !+ DIFF_V(:,:,I) |
1253 |
|
|
ETA_C1(:,:,I) = ETA_C2(:,:,I) !+ DIFF_ETA(:,:,I) |
1254 |
|
|
ENDDO |
1255 |
|
|
C---------------------------------------------------------- |
1256 |
|
|
ICONT = -1 |
1257 |
|
|
DO I = 1,nSxC |
1258 |
|
|
DO J = 1,nSyC |
1259 |
|
|
ICONT = ICONT + 1 |
1260 |
|
|
IndI = IM_C*(I-1) |
1261 |
|
|
IndJ = JM_C*(J-1) |
1262 |
|
|
|
1263 |
|
|
VAR_C1(:,:,:,:) = 0. |
1264 |
|
|
|
1265 |
|
|
J1 = 1+IndJ-OLY |
1266 |
|
|
J2 = JM_C+IndJ+OLY |
1267 |
|
|
|
1268 |
|
|
JJ1 = 1 |
1269 |
|
|
JJ2 = JM_C+OLY+OLY |
1270 |
|
|
|
1271 |
|
|
IF(1 +IndJ-OLY.LT.0) THEN |
1272 |
|
|
J1 = 1 |
1273 |
|
|
JJ1 = 4 |
1274 |
|
|
ENDIF |
1275 |
|
|
|
1276 |
|
|
IF(JM_C+IndJ+OLY.GT.NyC) THEN |
1277 |
|
|
J2 = NyC |
1278 |
|
|
JJ2 = JM_C |
1279 |
|
|
ENDIF |
1280 |
|
|
|
1281 |
|
|
VAR_C1(JJ1:JJ2,:,:,1) = U_C1(J1:J2,:,:) |
1282 |
|
|
VAR_C1(JJ1:JJ2,:,:,2) = V_C1(J1:J2,:,:) |
1283 |
|
|
VAR_C1(JJ1:JJ2,:,:,3) = T_C1(J1:J2,:,:) |
1284 |
|
|
VAR_C1(JJ1:JJ2,:,:,4) = S_C1(J1:J2,:,:) |
1285 |
|
|
VAR_C1(JJ1:JJ2,:,:,5) = ETA_C1(J1:J2,:,:) |
1286 |
|
|
|
1287 |
|
|
CALL MPI_SEND (VAR_C1, indexF, MPI_REAL8, |
1288 |
|
|
& MSTR_CHLD(NST_LEV)+ICONT, 3000, |
1289 |
|
|
& MPI_Comm_World,ierr) |
1290 |
|
|
|
1291 |
|
|
ENDDO |
1292 |
|
|
ENDDO |
1293 |
|
|
C--------------------------------------------------------------- |
1294 |
|
|
8888 CONTINUE |
1295 |
|
|
C-------------------------------------------------------- |
1296 |
|
|
C Close OBCs Files x CHILD MODEL |
1297 |
|
|
C-------------------------------------------------------- |
1298 |
|
|
C---------------------------------------------------- |
1299 |
|
|
C------------- MANDO SEGNALE DI OK AL CHILD |
1300 |
|
|
C---------------------------------------------------- |
1301 |
|
|
C-------------------------------------------------------- |
1302 |
|
|
C (1) READ FROM CHILD MODEL |
1303 |
|
|
C-------------------------------------------------------- |
1304 |
|
|
CALL MPI_RECV (TRANSPORT_WEST, 1, MPI_REAL8, |
1305 |
|
|
& MSTR_CHLD(NST_LEV), 3000, |
1306 |
|
|
& MPI_COMM_World, status,ierr) |
1307 |
|
|
|
1308 |
|
|
CALL MPI_RECV (TRANSPORT_EAST, 1, MPI_REAL8, |
1309 |
|
|
& MSTR_CHLD(NST_LEV), 3000, |
1310 |
|
|
& MPI_COMM_World, status,ierr) |
1311 |
|
|
C--------------------------------------------------------- |
1312 |
|
|
C--------------------------------------------------------- |
1313 |
|
|
|
1314 |
|
|
ICONT=1 |
1315 |
|
|
|
1316 |
|
|
DO WHILE(ICONT.le.nSxC*nSyC) |
1317 |
|
|
write(iUnit,*) |
1318 |
|
|
& 'CALL MPI_RECV 3-D var from CHILD, index3F=', index3F |
1319 |
|
|
from= MPI_ANY_SOURCE |
1320 |
|
|
CALL MPI_RECV (globalC3D_a,index3F, MPI_REAL8, |
1321 |
|
|
& from, 3000, MPI_COMM_World, status,ierr) |
1322 |
|
|
|
1323 |
|
|
ICONT=ICONT+1 |
1324 |
|
|
|
1325 |
|
|
whm=status(MPI_SOURCE)-MSTR_CHLD(NST_LEV)+1 |
1326 |
|
|
|
1327 |
|
|
CALL MPI_GET_COUNT(status,MPI_REAL8,st_count,ierr) |
1328 |
|
|
|
1329 |
|
|
globalC3D(1+IndI_C(whm):IM_C+IndI_C(whm), |
1330 |
|
|
& 1+IndJ_C(whm):JM_C+IndJ_C(whm),:,:)= |
1331 |
|
|
& globalC3D_a(:,:,:,:) |
1332 |
|
|
END DO |
1333 |
|
|
C----------------------------- |
1334 |
|
|
ICONT=1 |
1335 |
|
|
|
1336 |
|
|
DO WHILE(ICONT.le.nSxC*nSyC) |
1337 |
|
|
from= MPI_ANY_SOURCE |
1338 |
|
|
CALL MPI_RECV (globalC2D_a,index2F, MPI_REAL8, |
1339 |
|
|
& from, 3000, MPI_COMM_World, status,ierr) |
1340 |
|
|
|
1341 |
|
|
ICONT=ICONT+1 |
1342 |
|
|
|
1343 |
|
|
whm=status(MPI_SOURCE)-MSTR_CHLD(NST_LEV)+1 |
1344 |
|
|
|
1345 |
|
|
CALL MPI_GET_COUNT(status,MPI_REAL8,st_count,ierr) |
1346 |
|
|
|
1347 |
|
|
globalC2D(1+IndI_C(whm):IM_C+IndI_C(whm), |
1348 |
|
|
& 1+IndJ_C(whm):JM_C+IndJ_C(whm),:)= |
1349 |
|
|
& globalC2D_a(:,:,:) |
1350 |
|
|
END DO |
1351 |
|
|
C-- end if rank=0 |
1352 |
|
|
ENDIF |
1353 |
|
|
|
1354 |
|
|
CALL MPI_BCAST(globalC3D,NxC*NyC*NrC*n3dC,MPI_REAL8, |
1355 |
|
|
& 0,NEST_COMM,ierr) |
1356 |
|
|
CALL MPI_BCAST(globalC2D,NxC*NyC*4,MPI_REAL8, |
1357 |
|
|
& 0,NEST_COMM,ierr) |
1358 |
|
|
|
1359 |
|
|
2323 CONTINUE |
1360 |
|
|
C======================================================= |
1361 |
|
|
C (1) READ FROM CHILD MODEL |
1362 |
|
|
C======================================================= |
1363 |
|
|
|
1364 |
|
|
C======================================================= |
1365 |
|
|
C (2) INTERPOLATIONS |
1366 |
|
|
C======================================================= |
1367 |
|
|
|
1368 |
|
|
C 3D VAR |
1369 |
|
|
C-------- |
1370 |
|
|
DO iVar = 1,15 ! tipo di variabile |
1371 |
|
|
DO K = 1,NrP |
1372 |
|
|
DO J = 1,NyP |
1373 |
|
|
DO I = WesternB+1,EasternB-1 |
1374 |
|
|
VAR3D_P(I,J,K,iVar) = 0. ! inizializzo |
1375 |
|
|
C WesternB side |
1376 |
|
|
|
1377 |
|
|
AREA_VOL = 0. !can be area or volume depend on the variable |
1378 |
|
|
|
1379 |
|
|
SELECT CASE(iVar) |
1380 |
|
|
CASE(1,5,9) |
1381 |
|
|
I_START = 1 |
1382 |
|
|
I_END = 9 |
1383 |
|
|
I_STEP = 1 !3 |
1384 |
|
|
CASE(2,6,10) |
1385 |
|
|
I_START = 1 |
1386 |
|
|
I_END = 9 !3 |
1387 |
|
|
I_STEP = 1 |
1388 |
|
|
CASE DEFAULT |
1389 |
|
|
I_START = 1 |
1390 |
|
|
I_END = 9 |
1391 |
|
|
I_STEP = 1 |
1392 |
|
|
END SELECT |
1393 |
|
|
|
1394 |
|
|
DO II = I_START,I_END,I_STEP |
1395 |
|
|
|
1396 |
|
|
IF (I_C2P(II,I,J).eq.0) cycle |
1397 |
|
|
IF (J_C2P(II,I,J).eq.0) cycle |
1398 |
|
|
|
1399 |
|
|
Indx = I_C2P(II,I,J) |
1400 |
|
|
Jndx = J_C2P(II,I,J) |
1401 |
|
|
|
1402 |
|
|
SELECT CASE(iVar) |
1403 |
|
|
|
1404 |
|
|
CASE (1,5,9) |
1405 |
|
|
VAR3D_P(I,J,K,ivar) = VAR3D_P(I,J,K,iVar) + |
1406 |
|
|
& globalC3D(Indx,Jndx,K,iVar)* |
1407 |
|
|
$ RAW_C(Indx,Jndx)* |
1408 |
|
|
& hFacW_C(Indx,Jndx,K) |
1409 |
|
|
|
1410 |
|
|
CASE (2,6,10) |
1411 |
|
|
VAR3D_P(I,J,K,ivar) = VAR3D_P(I,J,K,iVar) + |
1412 |
|
|
& globalC3D(Indx,Jndx,K,iVar)* |
1413 |
|
|
$ RAS_C(Indx,Jndx)* |
1414 |
|
|
& hFacS_C(Indx,Jndx,K) |
1415 |
|
|
|
1416 |
|
|
CASE DEFAULT |
1417 |
|
|
VAR3D_P(I,J,K,ivar) = VAR3D_P(I,J,K,iVar) + |
1418 |
|
|
& globalC3D(Indx,Jndx,K,iVar)* |
1419 |
|
|
$ RAC_C(Indx,Jndx)* |
1420 |
|
|
& hFacC_C(Indx,Jndx,K) |
1421 |
|
|
|
1422 |
|
|
AREA_VOL = AREA_VOL + |
1423 |
|
|
& RAC_C(Indx,Jndx)* hFacC_C(Indx,Jndx,K) |
1424 |
|
|
|
1425 |
|
|
END SELECT |
1426 |
|
|
ENDDO |
1427 |
|
|
C----------------------------------------------- |
1428 |
|
|
C Make a volume average |
1429 |
|
|
C---------------------------------------------- |
1430 |
|
|
SELECT CASE(iVar) |
1431 |
|
|
CASE (1,5,9) |
1432 |
|
|
VAR3D_P(I,J,K,ivar) = |
1433 |
|
|
& VAR3D_P(I,J,K,iVar)* |
1434 |
|
|
& INV_VOL_W_P(I,J,K) |
1435 |
|
|
if (hFacW_P(I,J,K).eq.0.) VAR3D_P(I,J,K,ivar)=0. |
1436 |
|
|
|
1437 |
|
|
CASE (2,6,10) |
1438 |
|
|
VAR3D_P(I,J,K,ivar) = |
1439 |
|
|
& VAR3D_P(I,J,K,iVar)* |
1440 |
|
|
& INV_VOL_S_P(I,J,K) |
1441 |
|
|
if (hFacS_P(I,J,K).eq.0.) VAR3D_P(I,J,K,ivar)=0. |
1442 |
|
|
CASE DEFAULT |
1443 |
|
|
IF (AREA_VOL.ne.0.) THEN |
1444 |
|
|
VAR3D_P(I,J,K,ivar) = |
1445 |
|
|
& VAR3D_P(I,J,K,iVar)/AREA_VOL |
1446 |
|
|
ENDIF |
1447 |
|
|
if (hFacC_P(I,J,K).eq.0.) VAR3D_P(I,J,K,ivar)=0. |
1448 |
|
|
END SELECT |
1449 |
|
|
ENDDO |
1450 |
|
|
ENDDO |
1451 |
|
|
ENDDO |
1452 |
|
|
ENDDO |
1453 |
|
|
|
1454 |
|
|
C 2D VAR |
1455 |
|
|
C-------- |
1456 |
|
|
DO iVar = 1,4 |
1457 |
|
|
DO J = 1,NyP |
1458 |
|
|
DO I = WesternB+1,EasternB-1 |
1459 |
|
|
VAR2D_P(I,J,iVar) = 0. |
1460 |
|
|
AREA_VOL = 0. |
1461 |
|
|
DO II = 1,9 |
1462 |
|
|
IF (I_C2P(II,I,J).eq.0) cycle |
1463 |
|
|
IF (J_C2P(II,I,J).eq.0) cycle |
1464 |
|
|
|
1465 |
|
|
Indx = I_C2P(II,I,J) |
1466 |
|
|
Jndx = J_C2P(II,I,J) |
1467 |
|
|
|
1468 |
|
|
VAR2D_P(I,J,ivar) = VAR2D_P(I,J,iVar) + |
1469 |
|
|
& globalC2D(Indx,Jndx,iVar)* |
1470 |
|
|
$ RAC_C(Indx,Jndx)* |
1471 |
|
|
& DEEP_C(Indx,Jndx,1) |
1472 |
|
|
|
1473 |
|
|
|
1474 |
|
|
AREA_VOL = AREA_VOL + |
1475 |
|
|
& RAC_C(Indx,Jndx)* DEEP_C(Indx,Jndx,1) |
1476 |
|
|
|
1477 |
|
|
ENDDO |
1478 |
|
|
C----------------------------- |
1479 |
|
|
IF ((RAC_P(I,J)*DEEP_P(I,J,1)).ne.0.) then |
1480 |
|
|
c IF (AREA_VOL.ne.0.) then |
1481 |
|
|
|
1482 |
|
|
VAR2D_P(I,J,ivar) = |
1483 |
|
|
& VAR2D_P(I,J,iVar)/ |
1484 |
|
|
& RAC_P(I,J) |
1485 |
|
|
ENDIF |
1486 |
|
|
C---------------------------- |
1487 |
|
|
ENDDO |
1488 |
|
|
ENDDO |
1489 |
|
|
ENDDO |
1490 |
|
|
|
1491 |
|
|
IF(rank.EQ.0) THEN |
1492 |
|
|
C-------------------------------------------------------- |
1493 |
|
|
C Write Files for PARENT MODEL |
1494 |
|
|
C-------------------------------------------------------- |
1495 |
|
|
C write(iUnit,*) ' (*) Open Files for PARENT MODEL' |
1496 |
|
|
|
1497 |
|
|
7575 CONTINUE |
1498 |
|
|
C---------------------------------------------------- |
1499 |
|
|
C------------- MANDO SEGNALE DI OK AL PARENT |
1500 |
|
|
C---------------------------------------------------- |
1501 |
|
|
CALL MPI_SEND (TRANSPORT_WEST, 1, MPI_REAL8, |
1502 |
|
|
& MSTR_PRNT(NST_LEV), 3000, |
1503 |
|
|
& MPI_Comm_World,ierr) |
1504 |
|
|
|
1505 |
|
|
CALL MPI_SEND (TRANSPORT_EAST, 1, MPI_REAL8, |
1506 |
|
|
& MSTR_PRNT(NST_LEV), 3000, |
1507 |
|
|
& MPI_Comm_World,ierr) |
1508 |
|
|
|
1509 |
|
|
ENDIF |
1510 |
|
|
C--------------------------------------------------------- |
1511 |
|
|
VCONT=VCONTP(rank) |
1512 |
|
|
|
1513 |
|
|
DO I = vstart,vstop |
1514 |
|
|
DO J = 1,nSyP |
1515 |
|
|
VCONT = VCONT + 1 |
1516 |
|
|
IndI = IM_P*(I-1) |
1517 |
|
|
IndJ = JM_P*(J-1) |
1518 |
|
|
C----------------------------------------------------------- |
1519 |
|
|
DO iVar=1,15 |
1520 |
|
|
c CALL MPI_SEND (VAR3D_P(1+IndI:IM_P+IndI |
1521 |
|
|
c & ,1+IndJ:JM_P+IndJ,:,iVar), |
1522 |
|
|
c & index_var3D,MPI_REAL8,MSTR_PRNT(NST_LEV)+VCONT, |
1523 |
|
|
c & 3000,MPI_Comm_World,ierr) |
1524 |
|
|
localP3D_a(:,:,:) = |
1525 |
|
|
& VAR3D_P(1+IndI:IM_P+IndI,1+IndJ:JM_P+IndJ,:,iVar) |
1526 |
|
|
CALL MPI_SEND ( localP3D_a, index_var3D, MPI_REAL8, |
1527 |
|
|
& MSTR_PRNT(NST_LEV)+VCONT, |
1528 |
|
|
& 3000, MPI_Comm_World, ierr ) |
1529 |
|
|
|
1530 |
|
|
ENDDO |
1531 |
|
|
|
1532 |
|
|
DO iVar=1,4 |
1533 |
|
|
c CALL MPI_SEND (VAR2D_P(1+IndI:IM_P+IndI |
1534 |
|
|
c & ,1+IndJ:JM_P+IndJ,iVar), |
1535 |
|
|
c & index_var2D,MPI_REAL8,MSTR_PRNT(NST_LEV)+VCONT, |
1536 |
|
|
c & 3000,MPI_Comm_World,ierr) |
1537 |
|
|
localP2D_a(:,:) = |
1538 |
|
|
& VAR2D_P(1+IndI:IM_P+IndI,1+IndJ:JM_P+IndJ,iVar) |
1539 |
|
|
CALL MPI_SEND ( localP2D_a, index_var2D, MPI_REAL8, |
1540 |
|
|
& MSTR_PRNT(NST_LEV)+VCONT, |
1541 |
|
|
& 3000, MPI_Comm_World, ierr ) |
1542 |
|
|
ENDDO |
1543 |
|
|
|
1544 |
|
|
C----------------------------------------------------------- |
1545 |
|
|
ENDDO |
1546 |
|
|
ENDDO |
1547 |
|
|
|
1548 |
|
|
CALL MPI_BARRIER(NEST_COMM,ierr) |
1549 |
|
|
ENDDO |
1550 |
|
|
C--------------------------------------------------------- |
1551 |
|
|
C======================================================= |
1552 |
|
|
C END MAIN LOOP |
1553 |
|
|
C======================================================= |
1554 |
|
|
CLOSE( iUnit ) |
1555 |
|
|
CALL MPI_FINALIZE(ierr) |
1556 |
|
|
C--------------------------------------------------------- |
1557 |
|
|
|
1558 |
|
|
101 FORMAT (I1) |
1559 |
|
|
|
1560 |
|
|
STOP |
1561 |
|
|
END |