/[MITgcm]/MITgcm_contrib/nesting_sannino/nest_driver/main.F
ViewVC logotype

Contents of /MITgcm_contrib/nesting_sannino/nest_driver/main.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.2 - (show annotations) (download)
Mon Nov 29 14:24:14 2010 UTC (14 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +39 -7 lines
check for consistent number of nesting-steps between parent & child:
 if they agree, set number of nesting-steps in driver accordingly;
 if not, stop cleanly (call MPI_FINALIZE & stop).

1 C $Header: /u/gcmpack/MITgcm_contrib/nesting_sannino/nest_driver/main.F,v 1.1 2010/11/28 03:27:56 jmc Exp $
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 INTEGER nNestSteps, nNestStepsP, nNestStepsC
253 C----------------------------------------------------
254 C MPI starts here
255 C----------------------------------------------------
256 CALL MPI_Init(ierr)
257 CALL MPI_Comm_size(MPI_COMM_WORLD,size,ierr)
258 CALL MPI_Comm_rank(MPI_COMM_WORLD,rank,ierr)
259 npd=size-(NCPUs_PRNT+NCPUs_CHLD(1))
260 if(rank.lt.npd) color=0
261 CALL MPI_COMM_SPLIT (MPI_COMM_WORLD, color,0,
262 & NEST_COMM,ierr)
263
264 CALL MPI_Comm_size(NEST_COMM,isize,ierr)
265 CALL MPI_Comm_rank(NEST_COMM,irank,ierr)
266
267 C--------------------------------------------------------
268 C- change local dir to rank_N and open log file
269 CALL SETDIR( rank )
270 WRITE(fNam,'(A,I4.4)') 'STDlog.', rank
271 OPEN( iUnit, FILE=fNam, STATUS='unknown')
272 mLoop = 0
273 C--------------------------------------------------------
274 C COMPUTE MASTER VALUES
275 C--------------------------------------------------------
276 MSTR_DRV(1) = 0
277
278 MSTR_PRNT(1) = npd
279 MSTR_CHLD(1) = NCPUs_PRNT + npd
280
281 DO Count_Lev = 2, NST_LEV_TOT
282 MSTR_DRV(Count_Lev) = MSTR_CHLD(Count_Lev-1) +
283 & NCPUs_CHLD(Count_Lev - 1)
284
285 MSTR_CHLD(Count_Lev) = MSTR_DRV(Count_Lev) + 1
286 MSTR_PRNT(Count_Lev) = MSTR_CHLD(Count_Lev-1)
287 ENDDO
288
289 vstart = 1+rank*(nSxP/MSTR_PRNT(1))
290 vstop = (1+rank)*(nSxP/MSTR_PRNT(1))
291 VCONT = (nSxP/npd)*(nSyP)*rank-1
292 VCONTP(rank) = VCONT
293
294 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
295 C-- Print out nesting parameter:
296
297 WRITE(iUnit,'(A)') '// ==================================='
298 WRITE(iUnit,'(A)') '// NEST_DRIVER parameters :'
299 WRITE(iUnit,'(A)') '// ==================================='
300 WRITE(iUnit,*) 'NEST_DRIVER: rank =', rank, ' ; color =',color
301 WRITE(iUnit,*) 'NEST_DRIVER: size =', size, ' ; npd =', npd
302 WRITE(iUnit,*) 'NEST_DRIVER: irank =', irank, ' ; isize =', isize
303 WRITE(iUnit,*) 'NEST_DRIVER: vstart =', vstart
304 WRITE(iUnit,*) 'NEST_DRIVER: vstop =', vstop
305 WRITE(iUnit,*) 'NEST_DRIVER: VCONTP =', VCONTP(rank)
306
307 WRITE(iUnit,*) 'NEST_DRIVER: NST_LEV_TOT =', NST_LEV_TOT
308 WRITE(iUnit,*) 'NEST_DRIVER: NST_LEV =', NST_LEV
309 WRITE(iUnit,*) 'NEST_DRIVER: NCPUs_PRNT =', NCPUs_PRNT
310 WRITE(iUnit,*) 'NEST_DRIVER: NCPUs_CHLD =', NCPUs_CHLD
311 WRITE(iUnit,*) 'NEST_DRIVER: MSTR_DRV =', MSTR_DRV
312 WRITE(iUnit,*) 'NEST_DRIVER: MSTR_PRNT =', MSTR_PRNT
313 WRITE(iUnit,*) 'NEST_DRIVER: MSTR_CHLD =', MSTR_CHLD
314
315 C--------------------------------------------------------
316 C COMPUTE DOMAIN DECOMPOSITION
317 C--------------------------------------------------------
318 c if(rank.eq.0) then
319
320 IM_C = int(NxC/nSxC)
321 JM_C = int(NyC/nSyC)
322
323 IM_P = int(NxP/nSxP)
324 JM_P = int(NyP/nSyP)
325
326 ICONT = 0
327 DO I = 1,nSxP
328 DO J = 1,nSyP
329 ICONT = ICONT + 1
330 IndI_P(ICONT) = IM_P*(I-1)
331 IndJ_P(ICONT) = JM_P*(J-1)
332 END DO
333 END DO
334
335 ICONT = 0
336 DO I = 1,nSxC
337 DO J = 1,nSyC
338 ICONT = ICONT + 1
339 IndI_C(ICONT) = IM_C*(I-1)
340 IndJ_C(ICONT) = JM_C*(J-1)
341 END DO
342 END DO
343
344 index = 6*JM_P*NrP*5
345 index_var3D = IM_P*JM_P*NrP
346 index_var2D = IM_P*JM_P
347
348 indexF = (JM_C+OLY+OLY)*NrC*2*5
349 index3F = IM_C*JM_C*NrC*n3dC
350 index2F = IM_C*JM_C*4
351
352 ALLOCATE( globalPA(6,JM_P,NrP,5) )
353 ALLOCATE( localP3D_a(IM_P,JM_P,NrP) )
354 ALLOCATE( localP2D_a(IM_P,JM_P) )
355
356 ALLOCATE( VAR_C1((JM_C+OLY+OLY),NrC,2,5) )
357 ALLOCATE( globalC3D_a(IM_C,JM_C,NrC,n3dC) )
358 ALLOCATE( globalC2D_a(IM_C,JM_C,4))
359
360 IF ( rank.EQ.0 ) THEN
361 C--------------------------------------------------------
362 C WARNING
363 C--------------------------------------------------------
364 write(iUnit,*) '*************************************'
365 write(iUnit,*) ' have you update geometric files?'
366 write(iUnit,*) ' in ./CHILD e ./PARENT'
367 write(iUnit,*) '*************************************'
368 C--------------------------------------------------------
369 C PARENT MODEL
370 C--------------------------------------------------------
371 write(iUnit,*) ' [1] Read PARENT model geometry'
372 C----------------------------------------------------
373 C XC & YC
374 C----------------------------------------------------
375 MSIZE = NxP*NyP*WORDLENGTH
376
377 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
378 & file=trim(dirNEST)//'/PARENT/XC.data',
379 & form='unformatted')
380
381 read (1,REC=1) Xo_P(:,:)
382 close(1)
383
384 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
385 & file=trim(dirNEST)//'/PARENT/YC.data',
386 & form='unformatted')
387
388 read (1,REC=1) Yo_P(:,:)
389 close(1)
390 C----------------------------------------------------
391 C XG & YG
392 C----------------------------------------------------
393 MSIZE = NxP*NyP*WORDLENGTH
394
395 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
396 & file=trim(dirNEST)//'/PARENT/XG.data',
397 & form='unformatted')
398
399 read (1,REC=1) Xg_P(:,:)
400 close(1)
401
402 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
403 & file=trim(dirNEST)//'/PARENT/YG.data',
404 & form='unformatted')
405
406 read (1,REC=1) Yg_P(:,:)
407 close(1)
408 C----------------------------------------------------
409 C Yu
410 C----------------------------------------------------
411 DO J = 1,NyP
412 DO I = 1,NxP
413 c Xu_P(I,J) = Xg_P(I,J)
414 Yu_P(I,J) = Yo_P(I,J)
415 ENDDO
416 ENDDO
417 C----------------------------------------------------
418 C Xv & Yv
419 C----------------------------------------------------
420 DO J = 1,NyP
421 DO I = 1,NxP
422 Xv_P(I,J) = Xo_P(I,J)
423 Yv_P(I,J) = Yg_P(I,J)
424 ENDDO
425 ENDDO
426 C----------------------------------------------------
427 C hFacC
428 C----------------------------------------------------
429 MSIZE = NxP*NyP*NrP*WORDLENGTH
430
431 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
432 & file=trim(dirNEST)//'/PARENT/hFacC.data',
433 & form='unformatted')
434
435 read (1,REC=1) hFacC_P(:,:,:)
436 close(1)
437 C----------------------------------------------------
438 C hFacW
439 C----------------------------------------------------
440 MSIZE = NxP*NyP*NrP*WORDLENGTH
441
442 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
443 & file=trim(dirNEST)//'/PARENT/hFacW.data',
444 & form='unformatted')
445
446 read (1,REC=1) hFacW_P(:,:,:)
447 close(1)
448 C----------------------------------------------------
449 C hFacS
450 C----------------------------------------------------
451 MSIZE = NxP*NyP*NrP*WORDLENGTH
452
453 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
454 & file=trim(dirNEST)//'/PARENT/hFacS.data',
455 & form='unformatted')
456
457 read (1,REC=1) hFacS_P(:,:,:)
458 close(1)
459 C----------------------------------------------------
460 C RAC
461 C----------------------------------------------------
462 MSIZE = NxP*NyP*WORDLENGTH
463
464 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
465 & file=trim(dirNEST)//'/PARENT/RAC.data',
466 & form='unformatted')
467
468 read (1,REC=1) RAC_P(:,:)
469 close(1)
470 C----------------------------------------------------
471 C RAW
472 C----------------------------------------------------
473 MSIZE = NxP*NyP*WORDLENGTH
474
475 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
476 & file=trim(dirNEST)//'/PARENT/RAW.data',
477 & form='unformatted')
478
479 read (1,REC=1) RAW_P(:,:)
480 close(1)
481 C----------------------------------------------------
482 C RAS
483 C----------------------------------------------------
484 MSIZE = NxP*NyP*WORDLENGTH
485
486 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
487 & file=trim(dirNEST)//'/PARENT/RAS.data',
488 & form='unformatted')
489
490 read (1,REC=1) RAS_P(:,:)
491 close(1)
492 C----------------------------------------------------
493 C MASK x PARENT
494 C----------------------------------------------------
495 DO K = 1,NrP
496 DO J = 1,NyP
497 DO I = 1,NxP
498 DEEP_P(i,j,k) = 0.
499 IF (hFacC_P(i,j,k).ne.0) then
500 DEEP_P(I,J,K) = 1.
501 ENDIF
502 ENDDO
503 ENDDO
504 ENDDO
505 C----------------------------------------------------
506 C 1/Volume (C)
507 C----------------------------------------------------
508 DO K = 1,NrP
509 DO J = 1,NyP
510 DO I = 1,NxP
511 INV_VOL_C_P(I,J,K) = 1.
512 IF ((RAC_P(I,J)*hFacC_P(I,J,K)).ne.0.) THEN
513 INV_VOL_C_P(I,J,K) = 1./(RAC_P(I,J)*hFacC_P(I,J,K))
514 ENDIF
515 ENDDO
516 ENDDO
517 ENDDO
518 C----------------------------------------------------
519 C 1/Volume (W)
520 C----------------------------------------------------
521 DO K = 1,NrP
522 DO J = 1,NyP
523 DO I = 1,NxP
524 INV_VOL_W_P(I,J,K) = 1.
525 IF ((RAW_P(I,J)*hFacW_P(I,J,K)).ne.0.) THEN
526 INV_VOL_W_P(I,J,K) = 1./(RAW_P(I,J)*hFacW_P(I,J,K))
527 ENDIF
528 ENDDO
529 ENDDO
530 ENDDO
531 C----------------------------------------------------
532 C 1/Volume (S)
533 C----------------------------------------------------
534 DO K = 1,NrP
535 DO J = 1,NyP
536 DO I = 1,NxP
537 INV_VOL_S_P(I,J,K) = 1.
538 IF ((RAS_P(I,J)*hFacS_P(I,J,K)).ne.0.) THEN
539 INV_VOL_S_P(I,J,K) = 1./(RAS_P(I,J)*hFacS_P(I,J,K))
540 ENDIF
541 ENDDO
542 ENDDO
543 ENDDO
544 C--------------------------------------------------------
545 C CHILD MODEL
546 C--------------------------------------------------------
547 write(iUnit,*) ' [2] Read CHILD model geometry'
548 C----------------------------------------------------
549 C XC & YC
550 C----------------------------------------------------
551 MSIZE = NxC*NyC*WORDLENGTH
552
553 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
554 & file=trim(dirNEST)//'/CHILD/XC.data',
555 & form='unformatted')
556
557 read (1,REC=1) Xo_C(:,:)
558 close(1)
559
560 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
561 & file=trim(dirNEST)//'/CHILD/YC.data',
562 & form='unformatted')
563
564 read (1,REC=1) Yo_C(:,:)
565 close(1)
566 C----------------------------------------------------
567 C XG & YG
568 C----------------------------------------------------
569 MSIZE = NxC*NyC*WORDLENGTH
570
571 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
572 & file=trim(dirNEST)//'/CHILD/XG.data',
573 & form='unformatted')
574
575 read (1,REC=1) Xg_C(:,:)
576 close(1)
577
578 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
579 & file=trim(dirNEST)//'/CHILD/YG.data',
580 & form='unformatted')
581
582 read (1,REC=1) Yg_C(:,:)
583 close(1)
584 C----------------------------------------------------
585 C Xu & Yu
586 C----------------------------------------------------
587 DO J = 1,NyC
588 DO I = 1,NxC
589 Xu_C(I,J) = Xg_C(I,J)
590 Yu_C(I,J) = Yo_C(I,J)
591 ENDDO
592 ENDDO
593 C----------------------------------------------------
594 C Xv & Yv
595 C----------------------------------------------------
596 DO J = 1,NyC
597 DO I = 1,NxC
598 Xv_C(I,J) = Xo_C(I,J)
599 Yv_C(I,J) = Yg_C(I,J)
600 ENDDO
601 ENDDO
602 C----------------------------------------------------
603 C hFacC
604 C----------------------------------------------------
605 MSIZE = NxC*NyC*NrC*WORDLENGTH
606
607 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
608 & file=trim(dirNEST)//'/CHILD/hFacC.data',
609 & form='unformatted')
610
611 read (1,REC=1) hFacC_C(:,:,:)
612 close(1)
613 C----------------------------------------------------
614 C hFacW
615 C----------------------------------------------------
616 MSIZE = NxC*NyC*NrC*WORDLENGTH
617
618 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
619 & file=trim(dirNEST)//'/CHILD/hFacW.data',
620 & form='unformatted')
621
622 read (1,REC=1) hFacW_C(:,:,:)
623 close(1)
624 C----------------------------------------------------
625 C hFacC
626 C----------------------------------------------------
627 MSIZE = NxC*NyC*NrC*WORDLENGTH
628
629 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
630 & file=trim(dirNEST)//'/CHILD/hFacS.data',
631 & form='unformatted')
632
633 read (1,REC=1) hFacS_C(:,:,:)
634 close(1)
635 C----------------------------------------------------
636 C RAC
637 C----------------------------------------------------
638 MSIZE = NxC*NyC*WORDLENGTH
639
640 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
641 & file=trim(dirNEST)//'/CHILD/RAC.data',
642 & form='unformatted')
643
644 read (1,REC=1) RAC_C(:,:)
645 close(1)
646 C----------------------------------------------------
647 C RAW
648 C----------------------------------------------------
649 MSIZE = NxC*NyC*WORDLENGTH
650
651 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
652 & file=trim(dirNEST)//'/CHILD/RAW.data',
653 & form='unformatted')
654
655 read (1,REC=1) RAW_C(:,:)
656 close(1)
657 C----------------------------------------------------
658 C RAS
659 C----------------------------------------------------
660 MSIZE = NxC*NyC*WORDLENGTH
661
662 open(unit=1,ACCESS='direct',RECL=MSIZE,STATUS='OLD',
663 & file=trim(dirNEST)//'/CHILD/RAS.data',
664 & form='unformatted')
665
666 read (1,REC=1) RAS_C(:,:)
667 close(1)
668 C----------------------------------------------------
669 C MASK x CHILD
670 C----------------------------------------------------
671 DO K = 1,NrC
672 DO J = 1,NyC
673 DO I = 1,NxC
674 DEEP_C(i,j,k) = 0.
675 IF (hFacC_C(i,j,k).ne.0) then
676 DEEP_C(I,J,K) = 1.
677 ENDIF
678 ENDDO
679 ENDDO
680 ENDDO
681
682 C----------------------------------------------------
683 C __/________ ___________
684 C | | | | ||
685 C > o > | | |
686 C |__/__|_____|
687 C | | |
688 C > o > |
689 C |_____|_____|_____|_____|
690 C
691 C----------------------------------------------------
692 write(iUnit,*) ' [3] Compute J index P-->C'
693 C--------------------------------------------------------
694 C Compute J indicies for mapping P->C (I)
695 C--------------------------------------------------------
696 I = 1
697 II = WesternB
698
699 DO J = 1,NyC
700 P2C_U(J) = 0.
701 DO JJ = 1,NyP-1
702 YF = Yg_C(I,J)
703 YP1 = Yg_P(II,JJ)
704 YP2 = Yg_P(II,JJ+1)
705 IF (YF.ge.YP1.and.YF.lt.YP2) THEN
706 P2C_U(J) = JJ
707 ENDIF
708 ENDDO
709 ENDDO
710 C--------------------------------------------------------
711 C Compute J indicies for mapping P->C (II)
712 C--------------------------------------------------------
713 I = 1
714 II = WesternB
715
716 DO J = 1,NyC
717 P2C_linU(J) = 0.
718 DO JJ = 1,NyP-1
719 YF = Yu_C(I,J)
720 YP1 = Yu_P(II,JJ)
721 YP2 = Yu_P(II,JJ+1)
722 IF (YF.ge.YP1.and.YF.lt.YP2) THEN
723 P2C_linU(J) = JJ
724 ENDIF
725 ENDDO
726 ENDDO
727 C--------------------------------------------------------
728 C Compute J indicies for mapping P->C (III)
729 C--------------------------------------------------------
730 I = 1
731 II = WesternB
732
733 DO J = 1,NyC
734 DO JJ = 1,NyP-1
735 YF = Yu_C(I,J)
736 YP1 = Yu_P(II,JJ)
737 IF (YF.eq.YP1) THEN
738 WO3_linU(J) = 0
739 if (J+1.le.NyC) WO3_linU(J+1) = 1
740 if (J+2.le.NyC) WO3_linU(J+2) = 2
741 ENDIF
742 ENDDO
743 ENDDO
744 C--------------------Lower bound
745 DO J = 1,NyC
746 DO JJ = 1,NyP-1
747 YF = Yu_C(I,J)
748 YP1 = Yu_P(II,JJ)
749 IF (YF.eq.YP1) THEN
750 WO3_linU(J) = 0
751 if (J-1.gt.0) WO3_linU(J-1) = 2
752 if (J-2.gt.0) WO3_linU(J-2) = 1
753 GOTO 2345
754 ENDIF
755 ENDDO
756 ENDDO
757 2345 CONTINUE
758 C--------------------Upper bound
759 DO J = NyC,1,-1
760 DO JJ = 1,NyP-1
761 YF = Yu_C(I,J)
762 YP1 = Yu_P(II,JJ)
763 IF (YF.eq.YP1) THEN
764 WO3_linU(J) = 0
765 if (J+1.le.NyC) WO3_linU(J+1) = 1
766 if (J+2.le.NyC) WO3_linU(J+2) = 2
767 GOTO 2346
768 ENDIF
769 ENDDO
770 ENDDO
771 2346 CONTINUE
772 C--------------------------------------------------------
773 C Compute J indicies for mapping P->C (IV)
774 C--------------------------------------------------------
775 I = 1
776 II = WesternB
777
778 DO J = 1,NyC
779 P2C_linV(J) = 0.
780 DO JJ = 1,NyP-1
781 YF = Yv_C(I,J)
782 YP1 = Yv_P(II,JJ)
783 YP2 = Yv_P(II,JJ+1)
784 IF (YF.ge.YP1.and.YF.lt.YP2) THEN
785 P2C_linV(J) = JJ
786 ENDIF
787 ENDDO
788 ENDDO
789 C--------------------------------------------------------
790 C Compute J indicies for mapping P->C (V)
791 C--------------------------------------------------------
792 I = 1
793 II = WesternB
794
795 DO J = 1,NyC
796 DO JJ = 1,NyP-1
797 YF = Yv_C(I,J)
798 YP1 = Yv_P(II,JJ)
799 IF (YF.eq.YP1) THEN
800 WO3_linV(J) = 0
801 if (J+1.le.NyC) WO3_linV(J+1) = 1
802 if (J+2.le.NyC) WO3_linV(J+2) = 2
803 ENDIF
804 ENDDO
805 ENDDO
806 C--------------------Lower bound
807 DO J = 1,NyC
808 DO JJ = 1,NyP-1
809 YF = Yv_C(I,J)
810 YP1 = Yv_P(II,JJ)
811 IF (YF.eq.YP1) THEN
812 WO3_linV(J) = 0
813 if (J-1.gt.0) WO3_linV(J-1) = 2
814 if (J-2.gt.0) WO3_linV(J-2) = 1
815 GOTO 23451
816 ENDIF
817 ENDDO
818 ENDDO
819 23451 CONTINUE
820 C--------------------Upper bound
821 DO J = NyC,1,-1
822 DO JJ = 1,NyP-1
823 YF = Yv_C(I,J)
824 YP1 = Yv_P(II,JJ)
825 IF (YF.eq.YP1) THEN
826 WO3_linV(J) = 0
827 if (J+1.le.NyC) WO3_linV(J+1) = 1
828 if (J+2.le.NyC) WO3_linV(J+2) = 2
829 GOTO 23461
830 ENDIF
831 ENDDO
832 ENDDO
833 23461 CONTINUE
834 C--------------------------------------------------------
835 C Compute J indicies for mapping P->C (V)
836 C--------------------------------------------------------
837 write(iUnit,*) ' [5] Compute J index P-->C for (o)'
838 I = 1
839 II = WesternB
840
841 DO J = 1,NyC
842 P2C_o(J) = 0.
843 DO JJ = 1,NyP-1
844 YF = Yo_C(I,J)
845 YP1 = Yg_P(II,JJ)
846 YP2 = Yg_P(II,JJ+1)
847 IF (YF.gt.YP1.and.YF.lt.YP2) THEN
848 P2C_o(J) = JJ
849 ENDIF
850 ENDDO
851 ENDDO
852 C--------------------------------------------------------
853 C Compute J indicies for mapping P->C (VI)
854 C--------------------------------------------------------
855 write(iUnit,*) ' [6] Compute J index P-->C for (v bilinear)'
856 I = 1
857 II = WesternB
858
859 DO J = 1,NyC
860 DO JJ = 2,NyP-1
861 YF = Yv_C(I,J)
862 YP1 = Yv_P(II,JJ)
863 YP2 = Yv_P(II,JJ+1)
864 YP3 = Yv_P(II,JJ-1)
865
866 IF (YF.ge.YP1.and.YF.lt.YP2) THEN
867 P2C1_V(J) = JJ
868 P2C2_V(J) = JJ+1
869 ENDIF
870 ENDDO
871 ENDDO
872 C--------------------------------------------------------
873 C Look for the 9 CHILD indicies in PARENT grid cell
874 C--------------------------------------------------------
875 write(iUnit,*) ' [8] Compute I J index C-->P for (o)'
876
877 DO J = 1,NyP
878 DO I = 1,NxP
879 I_C2P(:,I,J) = 0
880 J_C2P(:,I,J) = 0
881
882 DO JJ = 1,NyC
883 DO II = 1,NxC
884 IF (Xo_C(II,JJ).eq.Xo_P(I,J).and.
885 & Yo_C(II,JJ).eq.Yo_P(I,J)) then
886
887 KK = 0
888 DO JJJ = JJ-1,JJ+1
889 DO III = II-1,II+1
890 KK = kk +1
891 if (III.lt.1.or.III.gt.NxC) cycle
892 if (JJJ.lt.1.or.JJJ.gt.NyC) cycle
893 I_C2P(KK,I,J) = III
894 J_C2P(KK,I,J) = JJJ
895 ENDDO
896 ENDDO
897 ENDIF
898
899 ENDDO
900 ENDDO
901
902 ENDDO
903 ENDDO
904 C-- end if rank=0
905 ENDIF
906 C--------------------------------------------------------
907 C Broadcast all the above variables
908 C--------------------------------------------------------
909 CALL MPI_BCAST(I_C2P,9*NxP*NyP,MPI_INTEGER,
910 & 0,NEST_COMM,ierr)
911 CALL MPI_BCAST(J_C2P,9*NxP*NyP,MPI_INTEGER,
912 & 0,NEST_COMM,ierr)
913
914 CALL MPI_BCAST(RAC_C,NxC*NyC,MPI_REAL,
915 & 0,NEST_COMM,ierr)
916 CALL MPI_BCAST(hFacC_C,NxC*NyC*NrC,MPI_REAL,
917 & 0,NEST_COMM,ierr)
918 CALL MPI_BCAST(INV_VOL_C_P,NxP*NyP*NrP,MPI_REAL,
919 & 0,NEST_COMM,ierr)
920
921 CALL MPI_BCAST(RAW_C,NxC*NyC,MPI_REAL,
922 & 0,NEST_COMM,ierr)
923 CALL MPI_BCAST(hFacW_C,NxC*NyC*NrC,MPI_REAL,
924 & 0,NEST_COMM,ierr)
925 CALL MPI_BCAST(INV_VOL_W_P,NxP*NyP*NrP,MPI_REAL,
926 & 0,NEST_COMM,ierr)
927
928 CALL MPI_BCAST(RAS_C,NxC*NyC,MPI_REAL,
929 & 0,NEST_COMM,ierr)
930 CALL MPI_BCAST(hFacS_C,NxC*NyC*NrC,MPI_REAL,
931 & 0,NEST_COMM,ierr)
932 CALL MPI_BCAST(INV_VOL_S_P,NxP*NyP*NrP,MPI_REAL,
933 & 0,NEST_COMM,ierr)
934
935 CALL MPI_BCAST(DEEP_C,NxC*NyC*NrC,MPI_REAL,
936 & 0,NEST_COMM,ierr)
937 CALL MPI_BCAST(RAC_P,NxP*NyP,MPI_REAL,
938 & 0,NEST_COMM,ierr)
939
940 CALL MPI_BCAST(IM_P,1,MPI_INTEGER,
941 & 0,NEST_COMM,ierr)
942 CALL MPI_BCAST(JM_P,1,MPI_INTEGER,
943 & 0,NEST_COMM,ierr)
944 CALL MPI_BCAST(index_var3D,1,MPI_INTEGER,
945 & 0,NEST_COMM,ierr)
946 CALL MPI_BCAST(index_var2D,1,MPI_INTEGER,
947 & 0,NEST_COMM,ierr)
948
949 CALL MPI_BCAST(DEEP_P,NxP*NyP*NrP,MPI_REAL,
950 & 0,NEST_COMM,ierr)
951 CALL MPI_BCAST(hFacS_P,NxP*NyP*NrP,MPI_REAL,
952 & 0,NEST_COMM,ierr)
953 CALL MPI_BCAST(hFacC_P,NxP*NyP*NrP,MPI_REAL,
954 & 0,NEST_COMM,ierr)
955 CALL MPI_BCAST(hFacW_P,NxP*NyP*NrP,MPI_REAL,
956 & 0,NEST_COMM,ierr)
957
958 C--------------------------------------------------------
959 if(rank.eq.0) then
960 C--------------------------------------------------------
961 DO K = 1,NrP
962 DO J = 1,NyP
963 DO I = WesternB+1,EasternB-1
964 C- WesternB side
965
966 DO II = 1,9
967 IF (I_C2P(II,I,J).eq.0) cycle
968 IF (J_C2P(II,I,J).eq.0) cycle
969
970 Indx = I_C2P(II,I,J)
971 Jndx = J_C2P(II,I,J)
972 ENDDO
973 ENDDO
974 ENDDO
975 ENDDO
976 C---------------------------------------------------------
977 ONOFF=0
978 endif
979 C---------------------------------------------------------
980 C Check parameter consistency across components:
981 C If inconsistent, send error code (-1) to every body and stop.
982 C- For now, just check number of nesting-steps between components
983 C---------------------------------------------------------
984 IF ( rank.EQ.0 ) THEN
985 C- Receive what the parent expects in term of nesting-exchanges Nb
986 CALL MPI_RECV( nNestStepsP, 1, MPI_INTEGER,
987 & MSTR_PRNT(NST_LEV), 3000,
988 & MPI_Comm_World, status, ierr )
989 C- Receive what the child expects in term of nesting-exchanges Nb
990 CALL MPI_RECV( nNestStepsC, 1, MPI_INTEGER,
991 & MSTR_CHLD(NST_LEV), 3000,
992 & MPI_Comm_World, status, ierr )
993 IF ( nNestStepsP .EQ. nNestStepsC ) THEN
994 nNestSteps = nNestStepsP
995 ELSE
996 WRITE(iUnit,'(A,I8)') ' ===== nNestStepsP =', nNestStepsP
997 WRITE(iUnit,'(A,I8)') ' ===== nNestStepsC =', nNestStepsC
998 nNestSteps = -1
999 ENDIF
1000 ENDIF
1001 c CALL MPI_BCAST( nNestSteps, 1, MPI_INTEGER,
1002 c & 0, NEST_COMM, ierr )
1003 C Note: above is redundant with following call
1004
1005 C- Broadcast error code (-1) from World-Master to every one in World
1006 C Note: better than from: MSTR_DRV(NST_LEV) since, to stop cleanly
1007 C with every one calling MPI_FINALIZE & stopping,
1008 C error needs to be sent to everybody.
1009 CALL MPI_BCAST( nNestSteps, 1, MPI_INTEGER,
1010 & 0, MPI_Comm_World, ierr )
1011 WRITE(iUnit,'(A,I8)') ' - - - nNestSteps =', nNestSteps
1012 C--------------------------------------------------------
1013 C BEGIN MAIN LOOP
1014 C--------------------------------------------------------
1015 DO mLoop=1,nNestSteps
1016 WRITE(iUnit,'(A,I8)') '== Main Loop , iter=', mLoop
1017 if(rank.eq.0) then
1018 C--------------------------------------------------------
1019 C (1) READ FROM PARENT MODEL
1020 C--------------------------------------------------------
1021 ICONT=1
1022 DO WHILE(ICONT.le.nSxP*nSyP)
1023 from= MPI_ANY_SOURCE
1024
1025 CALL MPI_RECV (globalPA, index, MPI_REAL8,
1026 & FROM, 3000,
1027 & MPI_COMM_World, status,ierr)
1028
1029 ICONT=ICONT+1
1030
1031 whm=status(MPI_SOURCE)-MSTR_PRNT(NST_LEV)+1
1032
1033 CALL MPI_GET_COUNT(status,MPI_REAL8,st_count,ierr)
1034
1035 DO II = 1,6
1036 IF (globalPA(II,1,1,1).ne.-999.) THEN
1037 globalP1(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1038 & globalPA(II,1:JM_P,:,1)
1039 globalP2(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1040 & globalPA(II,1:JM_P,:,2)
1041 globalP3(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1042 & globalPA(II,1:JM_P,:,3)
1043 globalP4(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1044 & globalPA(II,1:JM_P,:,4)
1045 globalP5(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1046 & globalPA(II,1:JM_P,:,5)
1047 ENDIF
1048
1049 ENDDO
1050 ENDDO
1051 C--------------------------------------------------------
1052 C Start interpolation for CHILD
1053 C--------------------------------------------------------
1054 CALL INTERPOLATION_P2C (
1055 & globalP1,globalP2,globalP3,globalP4,globalP5,
1056 & NxP,NyP,NrP,
1057 & NxC,NyC,NrC,
1058 $ WesternB,EasternB,
1059 $ P2C_U,P2C_V,P2C_o,P2C1_V,P2C2_V,P2C1_o,P2C2_o,
1060 $ P2C_linU,WO3_linU,P2C_linV,WO3_linV,
1061 $ Xv_C,Yv_C,Xv_P,Yv_P,
1062 $ T_C1,S_C1,U_C1,V_C1,ETA_C1,
1063 $ DEEP_C,DEEP_P
1064 & )
1065 C==============================================================
1066 C Open Files from PARENT MODEL
1067 C==============================================================
1068 ICONT=1
1069
1070 do while(ICONT.le.nSxP*nSyP)
1071 from= MPI_ANY_SOURCE
1072
1073 CALL MPI_RECV (globalPA, index, MPI_REAL8,
1074 & FROM, 3000,
1075 & MPI_COMM_World, status,ierr)
1076
1077 ICONT=ICONT+1
1078
1079 whm=status(MPI_SOURCE)-MSTR_PRNT(NST_LEV)+1
1080
1081 CALL MPI_GET_COUNT(status,MPI_REAL8,st_count,ierr)
1082
1083 DO II = 1,6
1084 IF (globalPA(II,1,1,1).ne.-999.) THEN
1085 globalP1(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1086 & globalPA(II,1:JM_P,:,1)
1087 globalP2(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1088 & globalPA(II,1:JM_P,:,2)
1089 globalP3(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1090 & globalPA(II,1:JM_P,:,3)
1091 globalP4(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1092 & globalPA(II,1:JM_P,:,4)
1093 globalP5(II,1+IndJ_P(whm):JM_P+IndJ_P(whm),:) =
1094 & globalPA(II,1:JM_P,:,5)
1095 ENDIF
1096 ENDDO
1097
1098 end do
1099 C--------------------------------------------------------
1100 C Start inteprolation for CHILD
1101 C--------------------------------------------------------
1102 CALL INTERPOLATION_P2C (
1103 & globalP1,globalP2,globalP3,globalP4,globalP5,
1104 & NxP,NyP,NrP,
1105 & NxC,NyC,NrC,
1106 $ WesternB,EasternB,
1107 $ P2C_U,P2C_V,P2C_o,P2C1_V,P2C2_V,P2C1_o,P2C2_o,
1108 $ P2C_linU,WO3_linU,P2C_linV,WO3_linV,
1109 $ Xv_C,Yv_C,Xv_P,Yv_P,
1110 $ T_C2,S_C2,U_C2,V_C2,ETA_C2,
1111 $ DEEP_C,DEEP_P
1112 & )
1113
1114 C==============================================================
1115 C Temporal Interpolation OBCs x CHILD MODEL
1116 C==============================================================
1117 C 0 1200
1118 C ---+--.--.--+---- Parent
1119 C
1120 C |--|--|--
1121 C 0 800
1122 C 400
1123 C------------------------------------------------------------
1124 DO I = 1,2 ! WesternB & EasternB
1125 DIFF_T(:,:,I) = (T_C2(:,:,I) - T_C1(:,:,I))/3.
1126 DIFF_S(:,:,I) = (S_C2(:,:,I) - S_C1(:,:,I))/3.
1127 DIFF_U(:,:,I) = (U_C2(:,:,I) - U_C1(:,:,I))/3.
1128 DIFF_V(:,:,I) = (V_C2(:,:,I) - V_C1(:,:,I))/3.
1129 DIFF_eta(:,:,I) = (eta_C2(:,:,I) - eta_C1(:,:,I))/3.
1130 ENDDO
1131 C-------------------------------------------------------------
1132 C Step 0 (Rec = 1 ==> WesternB)
1133 C------- (Rec = 2 ==> EasternB)
1134
1135 DO I = 1,2 !WesternB & EasternB
1136 T_C1(:,:,I) = T_C2(:,:,I) !+ DIFF_T(:,:,I)
1137 S_C1(:,:,I) = S_C2(:,:,I) !+ DIFF_S(:,:,I)
1138 U_C1(:,:,I) = U_C2(:,:,I) !+ DIFF_U(:,:,I)
1139 V_C1(:,:,I) = V_C2(:,:,I) !+ DIFF_V(:,:,I)
1140 ETA_C1(:,:,I) = ETA_C2(:,:,I) !+ DIFF_ETA(:,:,I)
1141 ENDDO
1142
1143 if(ONOFF.eq.0) then
1144 C---------------------------------------------------------------------
1145 ICONT = -1
1146 DO I = 1,nSxC
1147 DO J = 1,nSyC
1148 ICONT = ICONT + 1
1149 IndI = IM_C*(I-1)
1150 IndJ = JM_C*(J-1)
1151
1152 VAR_C1(:,:,:,:) = 0.
1153
1154 J1 = 1+IndJ-OLY
1155 J2 = JM_C+IndJ+OLY
1156
1157 JJ1 = 1
1158 JJ2 = JM_C+OLY+OLY
1159
1160 IF(1 +IndJ-OLY.LT.0) THEN
1161 J1 = 1
1162 JJ1 = 4
1163 ENDIF
1164
1165 IF(JM_C+IndJ+OLY.GT.NyC) THEN
1166 J2 = NyC
1167 JJ2 = JM_C
1168 ENDIF
1169
1170 VAR_C1(JJ1:JJ2,:,:,1) = U_C1(J1:J2,:,:)
1171 VAR_C1(JJ1:JJ2,:,:,2) = V_C1(J1:J2,:,:)
1172 VAR_C1(JJ1:JJ2,:,:,3) = T_C1(J1:J2,:,:)
1173 VAR_C1(JJ1:JJ2,:,:,4) = S_C1(J1:J2,:,:)
1174 VAR_C1(JJ1:JJ2,:,:,5) = ETA_C1(J1:J2,:,:)
1175
1176 CALL MPI_SEND (VAR_C1, indexF, MPI_REAL8,
1177 & MSTR_CHLD(NST_LEV)+ICONT, 3000,
1178 & MPI_Comm_World,ierr)
1179
1180 ENDDO
1181 ENDDO
1182 C----------------------------------------------------------------------
1183 c write(*,*) 'VIC: MANDO SEGNALE DI OK AL CHILD PER INIZIALIZZARE'
1184 ONOFF=1
1185 ENDIF
1186 c write(*,*) 'VIC: MANDO SEGNALE DI OK AL CHILD PER IL PASSO 1'
1187 C-----------------------------------------------------------------------
1188 ICONT = -1
1189 DO I = 1,nSxC
1190 DO J = 1,nSyC
1191 ICONT = ICONT + 1
1192 IndI = IM_C*(I-1)
1193 IndJ = JM_C*(J-1)
1194
1195 VAR_C1(:,:,:,:) = 0.
1196
1197 J1 = 1+IndJ-OLY
1198 J2 = JM_C+IndJ+OLY
1199
1200 JJ1 = 1
1201 JJ2 = JM_C+OLY+OLY
1202
1203 IF(1 +IndJ-OLY.LT.0) THEN
1204 J1 = 1
1205 JJ1 = 4
1206 ENDIF
1207
1208 IF(JM_C+IndJ+OLY.GT.NyC) THEN
1209 J2 = NyC
1210 JJ2 = JM_C
1211 ENDIF
1212
1213 VAR_C1(JJ1:JJ2,:,:,1) = U_C1(J1:J2,:,:)
1214 VAR_C1(JJ1:JJ2,:,:,2) = V_C1(J1:J2,:,:)
1215 VAR_C1(JJ1:JJ2,:,:,3) = T_C1(J1:J2,:,:)
1216 VAR_C1(JJ1:JJ2,:,:,4) = S_C1(J1:J2,:,:)
1217 VAR_C1(JJ1:JJ2,:,:,5) = ETA_C1(J1:J2,:,:)
1218
1219 CALL MPI_SEND (VAR_C1, indexF, MPI_REAL8,
1220 & MSTR_CHLD(NST_LEV)+ICONT, 3000,
1221 & MPI_Comm_World,ierr)
1222
1223 ENDDO
1224 ENDDO
1225 C---------------------------------------------------------------------
1226 goto 8888
1227
1228 C Step 1 (Rec = 3 ==> WesternB)
1229 C------- (Rec = 4 ==> EasternB)
1230
1231 DO I = 1,2 !WesternB & EasternB
1232 T_C1(:,:,I) = T_C2(:,:,I) !+ DIFF_T(:,:,I)
1233 S_C1(:,:,I) = S_C2(:,:,I) !+ DIFF_S(:,:,I)
1234 U_C1(:,:,I) = U_C2(:,:,I) !+ DIFF_U(:,:,I)
1235 V_C1(:,:,I) = V_C2(:,:,I) !+ DIFF_V(:,:,I)
1236 ETA_C1(:,:,I) = ETA_C2(:,:,I) !+ DIFF_ETA(:,:,I)
1237 ENDDO
1238 C----------------------------------------------------------
1239 ICONT = -1
1240 DO I = 1,nSxC
1241 DO J = 1,nSyC
1242 ICONT = ICONT + 1
1243 IndI = IM_C*(I-1)
1244 IndJ = JM_C*(J-1)
1245
1246 VAR_C1(:,:,:,:) = 0.
1247
1248 J1 = 1+IndJ-OLY
1249 J2 = JM_C+IndJ+OLY
1250
1251 JJ1 = 1
1252 JJ2 = JM_C+OLY+OLY
1253
1254 IF(1 +IndJ-OLY.LT.0) THEN
1255 J1 = 1
1256 JJ1 = 4
1257 ENDIF
1258
1259 IF(JM_C+IndJ+OLY.GT.NyC) THEN
1260 J2 = NyC
1261 JJ2 = JM_C
1262 ENDIF
1263
1264 VAR_C1(JJ1:JJ2,:,:,1) = U_C1(J1:J2,:,:)
1265 VAR_C1(JJ1:JJ2,:,:,2) = V_C1(J1:J2,:,:)
1266 VAR_C1(JJ1:JJ2,:,:,3) = T_C1(J1:J2,:,:)
1267 VAR_C1(JJ1:JJ2,:,:,4) = S_C1(J1:J2,:,:)
1268 VAR_C1(JJ1:JJ2,:,:,5) = ETA_C1(J1:J2,:,:)
1269
1270 CALL MPI_SEND (VAR_C1, indexF, MPI_REAL8,
1271 & MSTR_CHLD(NST_LEV)+ICONT, 3000,
1272 & MPI_Comm_World,ierr)
1273
1274 ENDDO
1275 ENDDO
1276 C-----------------------------------------------------------
1277 C Step 2 (Rec = 5 ==> WesternB)
1278 C------- (Rec = 6 ==> EasternB)
1279
1280 DO I = 1,2 !WesternB & EasternB
1281 T_C1(:,:,I) = T_C2(:,:,I) !+ DIFF_T(:,:,I)
1282 S_C1(:,:,I) = S_C2(:,:,I) !+ DIFF_S(:,:,I)
1283 U_C1(:,:,I) = U_C2(:,:,I) !+ DIFF_U(:,:,I)
1284 V_C1(:,:,I) = V_C2(:,:,I) !+ DIFF_V(:,:,I)
1285 ETA_C1(:,:,I) = ETA_C2(:,:,I) !+ DIFF_ETA(:,:,I)
1286 ENDDO
1287 C----------------------------------------------------------
1288 ICONT = -1
1289 DO I = 1,nSxC
1290 DO J = 1,nSyC
1291 ICONT = ICONT + 1
1292 IndI = IM_C*(I-1)
1293 IndJ = JM_C*(J-1)
1294
1295 VAR_C1(:,:,:,:) = 0.
1296
1297 J1 = 1+IndJ-OLY
1298 J2 = JM_C+IndJ+OLY
1299
1300 JJ1 = 1
1301 JJ2 = JM_C+OLY+OLY
1302
1303 IF(1 +IndJ-OLY.LT.0) THEN
1304 J1 = 1
1305 JJ1 = 4
1306 ENDIF
1307
1308 IF(JM_C+IndJ+OLY.GT.NyC) THEN
1309 J2 = NyC
1310 JJ2 = JM_C
1311 ENDIF
1312
1313 VAR_C1(JJ1:JJ2,:,:,1) = U_C1(J1:J2,:,:)
1314 VAR_C1(JJ1:JJ2,:,:,2) = V_C1(J1:J2,:,:)
1315 VAR_C1(JJ1:JJ2,:,:,3) = T_C1(J1:J2,:,:)
1316 VAR_C1(JJ1:JJ2,:,:,4) = S_C1(J1:J2,:,:)
1317 VAR_C1(JJ1:JJ2,:,:,5) = ETA_C1(J1:J2,:,:)
1318
1319 CALL MPI_SEND (VAR_C1, indexF, MPI_REAL8,
1320 & MSTR_CHLD(NST_LEV)+ICONT, 3000,
1321 & MPI_Comm_World,ierr)
1322
1323 ENDDO
1324 ENDDO
1325 C---------------------------------------------------------------
1326 8888 CONTINUE
1327 C--------------------------------------------------------
1328 C Close OBCs Files x CHILD MODEL
1329 C--------------------------------------------------------
1330 C----------------------------------------------------
1331 C------------- MANDO SEGNALE DI OK AL CHILD
1332 C----------------------------------------------------
1333 C--------------------------------------------------------
1334 C (1) READ FROM CHILD MODEL
1335 C--------------------------------------------------------
1336 CALL MPI_RECV (TRANSPORT_WEST, 1, MPI_REAL8,
1337 & MSTR_CHLD(NST_LEV), 3000,
1338 & MPI_COMM_World, status,ierr)
1339
1340 CALL MPI_RECV (TRANSPORT_EAST, 1, MPI_REAL8,
1341 & MSTR_CHLD(NST_LEV), 3000,
1342 & MPI_COMM_World, status,ierr)
1343 C---------------------------------------------------------
1344 C---------------------------------------------------------
1345
1346 ICONT=1
1347
1348 DO WHILE(ICONT.le.nSxC*nSyC)
1349 c write(iUnit,*)
1350 c & 'CALL MPI_RECV 3-D var from CHILD, index3F=', index3F
1351 from= MPI_ANY_SOURCE
1352 CALL MPI_RECV (globalC3D_a,index3F, MPI_REAL8,
1353 & from, 3000, MPI_COMM_World, status,ierr)
1354
1355 ICONT=ICONT+1
1356
1357 whm=status(MPI_SOURCE)-MSTR_CHLD(NST_LEV)+1
1358
1359 CALL MPI_GET_COUNT(status,MPI_REAL8,st_count,ierr)
1360
1361 globalC3D(1+IndI_C(whm):IM_C+IndI_C(whm),
1362 & 1+IndJ_C(whm):JM_C+IndJ_C(whm),:,:)=
1363 & globalC3D_a(:,:,:,:)
1364 END DO
1365 C-----------------------------
1366 ICONT=1
1367
1368 DO WHILE(ICONT.le.nSxC*nSyC)
1369 from= MPI_ANY_SOURCE
1370 CALL MPI_RECV (globalC2D_a,index2F, MPI_REAL8,
1371 & from, 3000, MPI_COMM_World, status,ierr)
1372
1373 ICONT=ICONT+1
1374
1375 whm=status(MPI_SOURCE)-MSTR_CHLD(NST_LEV)+1
1376
1377 CALL MPI_GET_COUNT(status,MPI_REAL8,st_count,ierr)
1378
1379 globalC2D(1+IndI_C(whm):IM_C+IndI_C(whm),
1380 & 1+IndJ_C(whm):JM_C+IndJ_C(whm),:)=
1381 & globalC2D_a(:,:,:)
1382 END DO
1383 C-- end if rank=0
1384 ENDIF
1385
1386 CALL MPI_BCAST(globalC3D,NxC*NyC*NrC*n3dC,MPI_REAL8,
1387 & 0,NEST_COMM,ierr)
1388 CALL MPI_BCAST(globalC2D,NxC*NyC*4,MPI_REAL8,
1389 & 0,NEST_COMM,ierr)
1390
1391 2323 CONTINUE
1392 C=======================================================
1393 C (1) READ FROM CHILD MODEL
1394 C=======================================================
1395
1396 C=======================================================
1397 C (2) INTERPOLATIONS
1398 C=======================================================
1399
1400 C 3D VAR
1401 C--------
1402 DO iVar = 1,15 ! tipo di variabile
1403 DO K = 1,NrP
1404 DO J = 1,NyP
1405 DO I = WesternB+1,EasternB-1
1406 VAR3D_P(I,J,K,iVar) = 0. ! inizializzo
1407 C WesternB side
1408
1409 AREA_VOL = 0. !can be area or volume depend on the variable
1410
1411 SELECT CASE(iVar)
1412 CASE(1,5,9)
1413 I_START = 1
1414 I_END = 9
1415 I_STEP = 1 !3
1416 CASE(2,6,10)
1417 I_START = 1
1418 I_END = 9 !3
1419 I_STEP = 1
1420 CASE DEFAULT
1421 I_START = 1
1422 I_END = 9
1423 I_STEP = 1
1424 END SELECT
1425
1426 DO II = I_START,I_END,I_STEP
1427
1428 IF (I_C2P(II,I,J).eq.0) cycle
1429 IF (J_C2P(II,I,J).eq.0) cycle
1430
1431 Indx = I_C2P(II,I,J)
1432 Jndx = J_C2P(II,I,J)
1433
1434 SELECT CASE(iVar)
1435
1436 CASE (1,5,9)
1437 VAR3D_P(I,J,K,ivar) = VAR3D_P(I,J,K,iVar) +
1438 & globalC3D(Indx,Jndx,K,iVar)*
1439 $ RAW_C(Indx,Jndx)*
1440 & hFacW_C(Indx,Jndx,K)
1441
1442 CASE (2,6,10)
1443 VAR3D_P(I,J,K,ivar) = VAR3D_P(I,J,K,iVar) +
1444 & globalC3D(Indx,Jndx,K,iVar)*
1445 $ RAS_C(Indx,Jndx)*
1446 & hFacS_C(Indx,Jndx,K)
1447
1448 CASE DEFAULT
1449 VAR3D_P(I,J,K,ivar) = VAR3D_P(I,J,K,iVar) +
1450 & globalC3D(Indx,Jndx,K,iVar)*
1451 $ RAC_C(Indx,Jndx)*
1452 & hFacC_C(Indx,Jndx,K)
1453
1454 AREA_VOL = AREA_VOL +
1455 & RAC_C(Indx,Jndx)* hFacC_C(Indx,Jndx,K)
1456
1457 END SELECT
1458 ENDDO
1459 C-----------------------------------------------
1460 C Make a volume average
1461 C----------------------------------------------
1462 SELECT CASE(iVar)
1463 CASE (1,5,9)
1464 VAR3D_P(I,J,K,ivar) =
1465 & VAR3D_P(I,J,K,iVar)*
1466 & INV_VOL_W_P(I,J,K)
1467 if (hFacW_P(I,J,K).eq.0.) VAR3D_P(I,J,K,ivar)=0.
1468
1469 CASE (2,6,10)
1470 VAR3D_P(I,J,K,ivar) =
1471 & VAR3D_P(I,J,K,iVar)*
1472 & INV_VOL_S_P(I,J,K)
1473 if (hFacS_P(I,J,K).eq.0.) VAR3D_P(I,J,K,ivar)=0.
1474 CASE DEFAULT
1475 IF (AREA_VOL.ne.0.) THEN
1476 VAR3D_P(I,J,K,ivar) =
1477 & VAR3D_P(I,J,K,iVar)/AREA_VOL
1478 ENDIF
1479 if (hFacC_P(I,J,K).eq.0.) VAR3D_P(I,J,K,ivar)=0.
1480 END SELECT
1481 ENDDO
1482 ENDDO
1483 ENDDO
1484 ENDDO
1485
1486 C 2D VAR
1487 C--------
1488 DO iVar = 1,4
1489 DO J = 1,NyP
1490 DO I = WesternB+1,EasternB-1
1491 VAR2D_P(I,J,iVar) = 0.
1492 AREA_VOL = 0.
1493 DO II = 1,9
1494 IF (I_C2P(II,I,J).eq.0) cycle
1495 IF (J_C2P(II,I,J).eq.0) cycle
1496
1497 Indx = I_C2P(II,I,J)
1498 Jndx = J_C2P(II,I,J)
1499
1500 VAR2D_P(I,J,ivar) = VAR2D_P(I,J,iVar) +
1501 & globalC2D(Indx,Jndx,iVar)*
1502 $ RAC_C(Indx,Jndx)*
1503 & DEEP_C(Indx,Jndx,1)
1504
1505
1506 AREA_VOL = AREA_VOL +
1507 & RAC_C(Indx,Jndx)* DEEP_C(Indx,Jndx,1)
1508
1509 ENDDO
1510 C-----------------------------
1511 IF ((RAC_P(I,J)*DEEP_P(I,J,1)).ne.0.) then
1512 c IF (AREA_VOL.ne.0.) then
1513
1514 VAR2D_P(I,J,ivar) =
1515 & VAR2D_P(I,J,iVar)/
1516 & RAC_P(I,J)
1517 ENDIF
1518 C----------------------------
1519 ENDDO
1520 ENDDO
1521 ENDDO
1522
1523 IF(rank.EQ.0) THEN
1524 C--------------------------------------------------------
1525 C Write Files for PARENT MODEL
1526 C--------------------------------------------------------
1527 C write(iUnit,*) ' (*) Open Files for PARENT MODEL'
1528
1529 7575 CONTINUE
1530 C----------------------------------------------------
1531 C------------- MANDO SEGNALE DI OK AL PARENT
1532 C----------------------------------------------------
1533 CALL MPI_SEND (TRANSPORT_WEST, 1, MPI_REAL8,
1534 & MSTR_PRNT(NST_LEV), 3000,
1535 & MPI_Comm_World,ierr)
1536
1537 CALL MPI_SEND (TRANSPORT_EAST, 1, MPI_REAL8,
1538 & MSTR_PRNT(NST_LEV), 3000,
1539 & MPI_Comm_World,ierr)
1540
1541 ENDIF
1542 C---------------------------------------------------------
1543 VCONT=VCONTP(rank)
1544
1545 DO I = vstart,vstop
1546 DO J = 1,nSyP
1547 VCONT = VCONT + 1
1548 IndI = IM_P*(I-1)
1549 IndJ = JM_P*(J-1)
1550 C-----------------------------------------------------------
1551 DO iVar=1,15
1552 c CALL MPI_SEND (VAR3D_P(1+IndI:IM_P+IndI
1553 c & ,1+IndJ:JM_P+IndJ,:,iVar),
1554 c & index_var3D,MPI_REAL8,MSTR_PRNT(NST_LEV)+VCONT,
1555 c & 3000,MPI_Comm_World,ierr)
1556 localP3D_a(:,:,:) =
1557 & VAR3D_P(1+IndI:IM_P+IndI,1+IndJ:JM_P+IndJ,:,iVar)
1558 CALL MPI_SEND ( localP3D_a, index_var3D, MPI_REAL8,
1559 & MSTR_PRNT(NST_LEV)+VCONT,
1560 & 3000, MPI_Comm_World, ierr )
1561
1562 ENDDO
1563
1564 DO iVar=1,4
1565 c CALL MPI_SEND (VAR2D_P(1+IndI:IM_P+IndI
1566 c & ,1+IndJ:JM_P+IndJ,iVar),
1567 c & index_var2D,MPI_REAL8,MSTR_PRNT(NST_LEV)+VCONT,
1568 c & 3000,MPI_Comm_World,ierr)
1569 localP2D_a(:,:) =
1570 & VAR2D_P(1+IndI:IM_P+IndI,1+IndJ:JM_P+IndJ,iVar)
1571 CALL MPI_SEND ( localP2D_a, index_var2D, MPI_REAL8,
1572 & MSTR_PRNT(NST_LEV)+VCONT,
1573 & 3000, MPI_Comm_World, ierr )
1574 ENDDO
1575
1576 C-----------------------------------------------------------
1577 ENDDO
1578 ENDDO
1579
1580 CALL MPI_BARRIER(NEST_COMM,ierr)
1581 ENDDO
1582 C---------------------------------------------------------
1583 C=======================================================
1584 C END MAIN LOOP
1585 C=======================================================
1586 CLOSE( iUnit )
1587 CALL MPI_FINALIZE(ierr)
1588 C---------------------------------------------------------
1589
1590 101 FORMAT (I1)
1591
1592 STOP
1593 END

  ViewVC Help
Powered by ViewVC 1.1.22