Newer
Older

WAUTELET Philippe
committed
!MNH_LIC Copyright 2005-2019 CNRS, Meteo-France and Universite Paul Sabatier
!MNH_LIC This is part of the Meso-NH software governed by the CeCILL-C licence

WAUTELET Philippe
committed
!MNH_LIC version 1. See LICENSE, CeCILL-C_V1-en.txt and CeCILL-C_V1-fr.txt
!MNH_LIC for details. version 1.
!-----------------------------------------------------------------
! ###############################
MODULE MODI_ADVEC_4TH_ORDER_AUX
! ###############################
!
INTERFACE
!
SUBROUTINE ADVEC_4TH_ORDER_ALGO(HLBCX, HLBCY, PFIELDT, KGRID, &
PMEANX, PMEANY,TPHALO2 )
!
USE MODE_ll
!
USE MODD_CONF
USE MODD_ARGSLIST_ll, ONLY : HALO2LIST_ll
!
CHARACTER (LEN=4), DIMENSION(2), INTENT(IN) :: HLBCX ! X direction LBC type
CHARACTER (LEN=4), DIMENSION(2), INTENT(IN) :: HLBCY ! Y direction LBC type
!

WAUTELET Philippe
committed
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PMEANX, PMEANY ! fluxes
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
REAL, DIMENSION(:,:,:), INTENT(IN) :: PFIELDT ! variable at t
INTEGER, INTENT(IN) :: KGRID ! C grid localisation
!
TYPE(HALO2_ll), OPTIONAL, POINTER :: TPHALO2 ! halo2 for the field at t
!
END SUBROUTINE ADVEC_4TH_ORDER_ALGO
!
!-------------------------------------------------------------------------------
!
FUNCTION MZF4(PA) RESULT(PMZF4)
!
REAL, DIMENSION(:,:,:), INTENT(IN) :: PA ! variable at flux
! side
REAL, DIMENSION(SIZE(PA,1),SIZE(PA,2),SIZE(PA,3)) :: PMZF4 ! result at mass
! localization
!
END FUNCTION MZF4
!
!-------------------------------------------------------------------------------
!
FUNCTION MZM4(PA) RESULT(PMZM4)
!
REAL, DIMENSION(:,:,:), INTENT(IN) :: PA ! variable at mass
! localization
REAL, DIMENSION(SIZE(PA,1),SIZE(PA,2),SIZE(PA,3)) :: PMZM4 ! result at flux
! localization
END FUNCTION MZM4
!
!-------------------------------------------------------------------------------
!
END INTERFACE
!
END MODULE MODI_ADVEC_4TH_ORDER_AUX
!
!-------------------------------------------------------------------------------
!
! ########################################################################
SUBROUTINE ADVEC_4TH_ORDER_ALGO(HLBCX, HLBCY, PFIELDT, KGRID, &
PMEANX, PMEANY,TPHALO2 )
! ########################################################################
!!
!!**** *ADVEC_4TH_ORDER_ALGO * - routine used to compute 4th order horizontal
!! advection fluxes of 3D prognostic variables
!!
!! PURPOSE
!! -------
!! The purpose of this routine is to compute 2sd or 4th order horizontal
!! advection fluxes of a prognostic variable.
!!
!!** METHOD
!! ------
!! In case of cyclic LBCs, the routine returns the scalar component of the
!! advection fluxes by applying a 4th order horizontal averaging operator to
!! the prognostic variable on each grid level. In the case of open LBCs, the
!! averaging operator degenerates to a 2nd order one on the first ring
!! inside the computationnal domain.
!! The "halo2" (or the second layer of the halo) of the prognostic
!! variable is passed as argument.
!!
!! IMPLICIT ARGUMENTS
!! ------------------
!!
!! MODULE MODD_ARGSLIST
!! HALO2LIST_ll : type for a list of "HALO2_lls"
!!
!! REFERENCE
!! ---------
!! Book2 of documentation ( routine ADVEC_4TH_ORDER )
!! User Interface for the MesoNH Parallel Package
!!
!! AUTHOR
!! ------
!! J.-P. Pinty * Laboratoire d'Aerologie*
!!
!! MODIFICATIONS
!! -------------
!! Original 25/10/05

ESCOBAR MUNOZ Juan
committed
!! Modif
!! J.Escobar 21/03/2013: for HALOK comment all NHALO=1 test
!!
!-------------------------------------------------------------------------------
!
!* 0. DECLARATIONS
! ------------
!

WAUTELET Philippe
committed
USE MODD_ARGSLIST_ll, ONLY: HALO2LIST_ll

WAUTELET Philippe
committed
!
USE MODE_ll
!
IMPLICIT NONE
!
!* 0.1 Declarations of dummy arguments :
!
CHARACTER (LEN=4), DIMENSION(2), INTENT(IN) :: HLBCX ! X direction LBC type
CHARACTER (LEN=4), DIMENSION(2), INTENT(IN) :: HLBCY ! Y direction LBC type
!

WAUTELET Philippe
committed
REAL, DIMENSION(:,:,:), INTENT(OUT) :: PMEANX, PMEANY ! fluxes
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
REAL, DIMENSION(:,:,:), INTENT(IN) :: PFIELDT ! variable at t
INTEGER, INTENT(IN) :: KGRID ! C grid localisation
!
TYPE(HALO2_ll), OPTIONAL, POINTER :: TPHALO2 ! halo2 for the field at t
!
!* 0.2 Declarations of local variables :
!
INTEGER:: IW,IE,IS,IN,IT,IB,IWF,IEF,ISF,INF ! Coordinate of forth order diffusion area
!
INTEGER:: IIB,IJB ! Begining useful area in x,y directions
INTEGER:: IIE,IJE ! End useful area in x,y directions
!
INTEGER:: ILUOUT,IRESP ! for prints
!
!-------------------------------------------------------------------------------
!
!* 0.3. COMPUTES THE DOMAIN DIMENSIONS
! ------------------------------
!
CALL GET_INDICE_ll(IIB,IJB,IIE,IJE)
!
!-------------------------------------------------------------------------------
!
!* 0.4. INITIALIZE THE FIELDS
! ---------------------
!
PMEANX(:,:,:) = 0.0
PMEANY(:,:,:) = 0.0
!
!-------------------------------------------------------------------------------
!
!
!* 1. CALCULATE THE NUMERICAL MEAN IN THE X DIRECTION
! -----------------------------------------------
!
SELECT CASE ( HLBCX(1) ) ! X direction LBC type: (1) for left side
!
!* 1.1 CYCLIC CASE IN THE X DIRECTION:
!
CASE ('CYCL') ! In that case one must have HLBCX(1) == HLBCX(2)
!
IW=IIB+1
IE=IIE
!
IF(KGRID == 2) THEN
IWF=IW-1
IEF=IE-1
ELSE
IWF=IW
IEF=IE
END IF
!
!* lateral boundary conditions
PMEANX(IWF-1,:,:) = (7.0*( PFIELDT(IW-1,:,:)+PFIELDT(IW-2,:,:) ) - &
( PFIELDT(IW,:,:)+TPHALO2%WEST(:,:) ) )/12.0
!
PMEANX(IEF+1,:,:) = (7.0*( PFIELDT(IE+1,:,:)+PFIELDT(IE,:,:) ) - &
( TPHALO2%EAST(:,:)+PFIELDT(IE-1,:,:) ) )/12.0
!
!* inner domain
PMEANX(IWF:IEF,:,:) = (7.0*( PFIELDT(IW:IE,:,:)+PFIELDT(IW-1:IE-1,:,:) ) - &
( PFIELDT(IW+1:IE+1,:,:)+PFIELDT(IW-2:IE-2,:,:) ) )/12.0
!
!!$!
!!$
!!$ IF(NHALO == 1) THEN
!!$ PMEANX(IWF-1,:,:) = (7.0*( PFIELDT(IW-1,:,:)+PFIELDT(IW-2,:,:) ) - &
!!$ ( PFIELDT(IW,:,:)+TPHALO2%WEST(:,:) ) )/12.0
!!$!
!!$ PMEANX(IEF+1,:,:) = (7.0*( PFIELDT(IE+1,:,:)+PFIELDT(IE,:,:) ) - &
!!$ ( TPHALO2%EAST(:,:)+PFIELDT(IE-1,:,:) ) )/12.0
!!$ ENDIF
!!$!
!!$ PMEANX(IWF:IEF,:,:) = (7.0*( PFIELDT(IW:IE,:,:)+PFIELDT(IW-1:IE-1,:,:) ) - &
!!$ ( PFIELDT(IW+1:IE+1,:,:)+PFIELDT(IW-2:IE-2,:,:) ) )/12.0
!!$!
!* 1.2 NON CYCLIC CASE IN THE X DIRECTION
!
CASE ('OPEN','WALL','NEST')
!
IF (LWEST_ll()) THEN
IF(KGRID == 2) THEN
IW=IIB+2 ! special case of C grid
ELSE
IW=IIB+1
END IF
ELSE

ESCOBAR MUNOZ Juan
committed
!!$ IF(NHALO == 1) THEN

ESCOBAR MUNOZ Juan
committed
!!$ ELSE
!!$ IW=IIB
!!$ ENDIF

ESCOBAR MUNOZ Juan
committed
!!$ IF (LEAST_ll() .OR. NHALO == 1) THEN
IF (LEAST_ll() ) THEN
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
! T. Maric
! IE=IIE-1 ! original
IE=IIE
ELSE
IE=IIE
END IF
!
IF(KGRID == 2) THEN
IWF=IW-1
IEF=IE-1
ELSE
IWF=IW
IEF=IE
END IF
!
! T. Maric. 16.1.2006.
! write(*,*)' IW, IE, IWF, IEF = ',IW, IE, IWF, IEF
! stop 'Stopping in advec_4th_order_aux.f90'
!
!* Use a second order scheme at the physical border
!
IF (LWEST_ll()) THEN
PMEANX(IWF-1,:,:) = 0.5*( PFIELDT(IW-1,:,:)+PFIELDT(IW-2,:,:) )
! T. Maric
! PMEANX(1,:,:) = PMEANX(IWF-1,:,:)
! extrapolate
!PMEANX(1,:,:) = 0.5*(3.0*PFIELDT(1,:,:) - PFIELDT(2,:,:))

ESCOBAR MUNOZ Juan
committed
!!$ ELSE IF (NHALO == 1) THEN
ELSE
PMEANX(IWF-1,:,:) = (7.0*( PFIELDT(IW-1,:,:)+PFIELDT(IW-2,:,:) ) - &
( PFIELDT(IW,:,:)+TPHALO2%WEST(:,:) ) )/12.0
ENDIF
!
IF (LEAST_ll()) THEN
PMEANX(IEF+1,:,:) = 0.5*( PFIELDT(IE+1,:,:)+PFIELDT(IE,:,:) )

ESCOBAR MUNOZ Juan
committed
!!$ ELSEIF (NHALO == 1) THEN
ELSE
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
PMEANX(IEF+1,:,:) = (7.0*( PFIELDT(IE+1,:,:)+PFIELDT(IE,:,:) ) - &
( TPHALO2%EAST(:,:)+PFIELDT(IE-1,:,:) ) )/12.0
ENDIF
!
!* Use a fourth order scheme elsewhere
!
PMEANX(IWF:IEF,:,:) = (7.0*( PFIELDT(IW:IE,:,:)+PFIELDT(IW-1:IE-1,:,:) ) - &
( PFIELDT(IW+1:IE+1,:,:)+PFIELDT(IW-2:IE-2,:,:) ) )/12.0
END SELECT
!
!-------------------------------------------------------------------------------
!
!* 2. COMPUTES THE 4TH ORDER MEAN IN THE Y DIRECTION
! ----------------------------------------------
!
IF ( .NOT. L2D ) THEN
SELECT CASE ( HLBCY(1) ) ! Y direction LBC type: (1) for left side
!
!* 2.1 CYCLIC CASE IN THE Y DIRECTION:
!
CASE ('CYCL') ! In that case one must have HLBCY(1) == HLBCY(2)
!
!
IS=IJB+1
IN=IJE
!
IF(KGRID == 3) THEN
ISF=IS-1
INF=IN-1
ELSE
ISF=IS
INF=IN
END IF
!
!* lateral boundary conditions
PMEANY(:,ISF-1,:) = (7.0*( PFIELDT(:,IS-1,:)+PFIELDT(:,IS-2,:) ) - &
( PFIELDT(:,IS,:)+TPHALO2%SOUTH(:,:) ) )/12.0
!
PMEANY(:,INF+1,:) = (7.0*( PFIELDT(:,IN+1,:)+PFIELDT(:,IN,:) ) - &
( TPHALO2%NORTH(:,:)+PFIELDT(:,IN-1,:) ) )/12.0
!
!* inner domain
PMEANY(:,ISF:INF,:) = (7.0*( PFIELDT(:,IS:IN,:)+PFIELDT(:,IS-1:IN-1,:)) - &
( PFIELDT(:,IS+1:IN+1,:)+PFIELDT(:,IS-2:IN-2,:) ))/12.0
!!$!
!!$ IF(NHALO == 1) THEN
!!$ PMEANY(:,ISF-1,:) = (7.0*( PFIELDT(:,IS,:)+PFIELDT(:,IS-1,:) ) - &
!!$ ( PFIELDT(:,IS+1,:)+TPHALO2%SOUTH(:,:) ) )/12.0
!!$!
!!$ PMEANY(:,ISF+1,:) = (7.0*( PFIELDT(:,IS,:)+PFIELDT(:,IS-1,:) ) - &
!!$ ( TPHALO2%NORTH(:,:)+PFIELDT(:,IS-2,:) ) )/12.0
!!$ ENDIF
!!$!
!!$ PMEANY(:,ISF:INF,:) = (7.0*( PFIELDT(:,IS:IN,:)+PFIELDT(:,IS-1:IN-1,:)) - &
!!$ ( PFIELDT(:,IS+1:IN+1,:)+PFIELDT(:,IS-2:IN-2,:) ))/12.0
!!$!
!* 2.2 NON CYCLIC CASE IN THE Y DIRECTION
!
CASE ('OPEN','WALL','NEST')
!
IF (LSOUTH_ll()) THEN
IF(KGRID == 3) THEN
IS=IJB+2 ! special case of C grid
ELSE
IS=IJB+1
END IF
ELSE

ESCOBAR MUNOZ Juan
committed
!!$ IF(NHALO == 1) THEN

ESCOBAR MUNOZ Juan
committed
!!$ ELSE
!!$ IS=IJB
!!$ ENDIF

ESCOBAR MUNOZ Juan
committed
!!$ IF (LNORTH_ll() .OR. NHALO == 1) THEN
IF (LNORTH_ll()) THEN
! T. Maric
! IN=IJE-1 ! original
IN=IJE
ELSE
IN=IJE
END IF
!
IF(KGRID == 3) THEN
ISF=IS-1
INF=IN-1
ELSE
ISF=IS
INF=IN
END IF
!
!* Use a second order scheme at the physical border
!
IF (LSOUTH_ll()) THEN
PMEANY(:,ISF-1,:) = 0.5*( PFIELDT(:,IS-1,:)+PFIELDT(:,IS-2,:) )
! T. Maric
! PMEANY(:,1,:) = PMEANY(:,ISF-1,:)
! extrapolate
!PMEANY(:,1,:) = 0.5*(3.0*PFIELDT(:,1,:) - PFIELDT(:,2,:))

ESCOBAR MUNOZ Juan
committed
!!$ ELSEIF (NHALO == 1) THEN
ELSE
!!$ PMEANY(:,ISF-1,:) = (7.0*( PFIELDT(:,IS,:)+PFIELDT(:,IS-1,:)) - &
!!$ ( PFIELDT(:,IS+1,:)+TPHALO2%SOUTH(:,:) ))/12.0
PMEANY(:,ISF-1,:) = (7.0*( PFIELDT(:,IS-1,:)+PFIELDT(:,IS-2,:)) - &
( PFIELDT(:,IS,:)+TPHALO2%SOUTH(:,:) ))/12.0
ENDIF
!
IF (LNORTH_ll()) THEN
PMEANY(:,INF+1,:) = 0.5*( PFIELDT(:,IN+1,:)+PFIELDT(:,IN,:) )

ESCOBAR MUNOZ Juan
committed
!!$ ELSEIF (NHALO == 1) THEN
ELSE
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
!!$ PMEANY(:,INF+1,:) = (7.0*( PFIELDT(:,IN,:)+PFIELDT(:,IN-1,:)) - &
!!$ ( TPHALO2%NORTH(:,:)+PFIELDT(:,IN-2,:) ))/12.0
PMEANY(:,INF+1,:) = (7.0*( PFIELDT(:,IN+1,:)+PFIELDT(:,IN,:)) - &
( TPHALO2%NORTH(:,:)+PFIELDT(:,IN-1,:) ))/12.0
ENDIF
!
!* Use a fourth order scheme elsewhere
!
PMEANY(:,ISF:INF,:) = (7.0*( PFIELDT(:,IS:IN,:)+PFIELDT(:,IS-1:IN-1,:)) - &
( PFIELDT(:,IS+1:IN+1,:)+PFIELDT(:,IS-2:IN-2,:) ))/12.0
!
END SELECT
ELSE
PMEANY(:,:,:) = 0.0
ENDIF
!
!-------------------------------------------------------------------------------
!
END SUBROUTINE ADVEC_4TH_ORDER_ALGO
!
!-------------------------------------------------------------------------------
!
! ################################
FUNCTION MZF4(PA) RESULT(PMZF4)
! ################################
!
!!**** *MZF4* - 4th order Shuman operator : mean operator in z direction for a
!! variable at a flux side
!!
!! PURPOSE
!! -------
!! The purpose of this function is to compute a 4th order mean value
!! along the z direction (K index) for a field PA localized at a z-flux
!! point (w point). The result is localized at a mass point.
!
!!** METHOD
!! ------
!! The result PMZF4(:,:,k) is defined by
!! PMZF4(:,:,k)=0.5*(PA(:,:,k)+PA(:,:,k+1)) at k=1 and size(PA,3)-1
!! PMZF4(:,:,k)=-999. at k=size(PA,3)
!! PMZF4(:,:,k)=7/12*(PA(:,:,k)+PA(:,:,k+1))
!! -1/12*(PA(:,:,k-1)+PA(:,:,k+2)) elsewhere
!!
!! EXTERNAL
!! --------
!! NONE
!!
!! IMPLICIT ARGUMENTS
!! ------------------
!! NONE
!!
!! REFERENCE
!! ---------
!! Book2 of documentation of Meso-NH (SHUMAN operators)
!! Technical specifications Report of The Meso-NH (chapters 3)
!!
!! AUTHOR
!! ------
!! J.-P. Pinty * Lab Aerologie *
!!
!! MODIFICATIONS
!! -------------
!! Original 25/10/05
!!
!-------------------------------------------------------------------------------
!
!* 0. DECLARATIONS
! ------------
!
IMPLICIT NONE
!
!* 0.1 Declarations of argument and result
!
!
REAL, DIMENSION(:,:,:), INTENT(IN) :: PA ! variable at flux
! side
REAL, DIMENSION(SIZE(PA,1),SIZE(PA,2),SIZE(PA,3)) :: PMZF4 ! result at mass
! localization
!
!* 0.2 Declarations of local variables
!
!
INTEGER :: JK ! loop index in z direction
INTEGER :: IKU ! upper bound in z direction of PA
!
INTEGER :: IIU,IJU,IIJU ! upper bounds in the x and y directions of PA
INTEGER :: JIJ,JIJK ! running loop indexes after linearisation
INTEGER :: JIJKOR1,JIJKEND1 ! loop boundaries
INTEGER :: JIJKOR2,JIJKEND2 ! loop boundaries
INTEGER :: JIJKOR3,JIJKEND3 ! loop boundaries
!
!-------------------------------------------------------------------------------
!
!* 1. DEFINITION OF MZF4
! ------------------
!
IIU = SIZE(PA,1)
IJU = SIZE(PA,2)
IKU = SIZE(PA,3)
!
IIJU = IIU*IJU
!
JIJKOR1 = 1 + IIJU
JIJKEND1 = 2*IIJU
!
!CDIR NODEP
!OCL NOVREC
DO JIJK=JIJKOR1 , JIJKEND1
PMZF4(JIJK-IIJU,1,1) = 0.5*( PA(JIJK-IIJU,1,1)+PA(JIJK,1,1) )
END DO
!
JIJKOR2 = 1 + JIJKEND1
JIJKEND2 = IIJU*IKU - IIJU
!
!CDIR NODEP
!OCL NOVREC
DO JIJK=JIJKOR2 , JIJKEND2
PMZF4(JIJK-IIJU,1,1) = (7.0*( PA(JIJK,1,1)+PA(JIJK-IIJU,1,1) ) - &
( PA(JIJK+IIJU,1,1)+PA(JIJK-2*IIJU,1,1) ) )/12.0
END DO
!
JIJKOR3 = 1 + JIJKEND2
JIJKEND3 = IIJU*IKU
!
!CDIR NODEP
!OCL NOVREC
DO JIJK=JIJKOR3 , JIJKEND3
PMZF4(JIJK-IIJU,1,1) = 0.5*( PA(JIJK-IIJU,1,1)+PA(JIJK,1,1) )
END DO
!
!CDIR NODEP
!OCL NOVREC
DO JIJ=1,IIJU
PMZF4(JIJ,1,IKU) = -999.
END DO
!
!-------------------------------------------------------------------------------
!
END FUNCTION MZF4
!
!-------------------------------------------------------------------------------
!
! ################################
FUNCTION MZM4(PA) RESULT(PMZM4)
! ################################
!
!!**** *MZM4* - 4th order Shuman operator : mean operator in z direction for a
!! mass variable
!!
!! PURPOSE
!! -------
!! The purpose of this function is to compute a 4th order mean value
!! along the z direction (K index) for a field PA localized at a mass
!! point. The result is localized at a z-flux point (w point).
!!
!!** METHOD
!! ------
!! The result PMZM4(:,:,k) is defined by
!! PMZM4(:,:,k)=0.5*(PA(:,:,k)+PA(:,:,k+1)) at k=2 and size(PA,3)
!! PMZM4(:,:,k)=-999. at k=1
!! PMZM4(:,:,k)=7/12*(PA(:,:,k)+PA(:,:,k+1))
!! -1/12*(PA(:,:,k-1)+PA(:,:,k+2)) elsewhere
!!
!! EXTERNAL
!! --------
!! NONE
!!
!! IMPLICIT ARGUMENTS
!! ------------------
!! NONE
!!
!! REFERENCE
!! ---------
!! Book2 of documentation of Meso-NH (SHUMAN operators)
!! Technical specifications Report of The Meso-NH (chapters 3)
!!
!! AUTHOR
!! ------
!! J.-P. Pinty * Lab Aerologie *
!!
!! MODIFICATIONS
!! -------------
!! Original 25/10/05
!!
!-------------------------------------------------------------------------------
!
!* 0. DECLARATIONS
! ------------
!
IMPLICIT NONE
!
!* 0.1 Declarations of argument and result
!
!
REAL, DIMENSION(:,:,:), INTENT(IN) :: PA ! variable at mass
! localization
REAL, DIMENSION(SIZE(PA,1),SIZE(PA,2),SIZE(PA,3)) :: PMZM4 ! result at flux
! localization
!
!* 0.2 Declarations of local variables
!
!
INTEGER :: JK ! loop index in z direction
INTEGER :: IKU ! upper bound in z direction of PA
!
INTEGER :: IIU,IJU,IIJU ! upper bounds in the x and y directions of PA
INTEGER :: JIJ,JIJK ! running loop indexes after linearisation
INTEGER :: JIJKOR1,JIJKEND1 ! loop boundaries
INTEGER :: JIJKOR2,JIJKEND2 ! loop boundaries
!
!-------------------------------------------------------------------------------
!
!* 1. DEFINITION OF MZM4
! ------------------
!
IIU = SIZE(PA,1)
IJU = SIZE(PA,2)
IKU = SIZE(PA,3)
!
IIJU = IIU*IJU
!
JIJKOR1 = 1 + IIJU
JIJKEND1 = JIJKOR1 + IIJU
!
!CDIR NODEP
!OCL NOVREC
DO JIJK=JIJKOR1 , JIJKEND1
PMZM4(JIJK,1,1) = 0.5*( PA(JIJK,1,1)+PA(JIJK-IIJU,1,1) )
END DO
!
JIJKOR2 = 1 + JIJKEND1
JIJKEND2 = IIJU*IKU - IIJU
!
!CDIR NODEP
!OCL NOVREC
DO JIJK=JIJKOR2 , JIJKEND2
PMZM4(JIJK,1,1) = (7.0*( PA(JIJK,1,1)+PA(JIJK-IIJU,1,1) ) - &
( PA(JIJK+IIJU,1,1)+PA(JIJK-2*IIJU,1,1) ) )/12.0
END DO
!
!CDIR NODEP
!OCL NOVREC
DO JIJ=1,IIJU
PMZM4(JIJ,1,IKU) = 0.5*( PA(JIJ,1,IKU)+PA(JIJ-IIJU,1,IKU) )
END DO
!
!CDIR NODEP
!OCL NOVREC
DO JIJ=1,IIJU
PMZM4(JIJ,1,1) = -999.
END DO
!
!-------------------------------------------------------------------------------
!
END FUNCTION MZM4