ROMS
Loading...
Searching...
No Matches
ana_m2clima.h
Go to the documentation of this file.
1!!
2 SUBROUTINE ana_m2clima (ng, tile, model)
3!
4!! git $Id$
5!!======================================================================
6!! Copyright (c) 2002-2025 The ROMS Group !
7!! Licensed under a MIT/X style license !
8!! See License_ROMS.md !
9!=======================================================================
10! !
11! This routine sets analytical 2D momentum climatology fields. !
12! !
13!=======================================================================
14!
15 USE mod_param
16 USE mod_ncparam
17!
18! Imported variable declarations.
19!
20 integer, intent(in) :: ng, tile, model
21!
22! Local variable declarations.
23!
24 character (len=*), parameter :: MyFile = &
25 & __FILE__
26!
27#include "tile.h"
28!
29 CALL ana_m2clima_tile (ng, tile, model, &
30 & lbi, ubi, lbj, ubj, &
31 & imins, imaxs, jmins, jmaxs)
32!
33! Set analytical header file name used.
34!
35#ifdef DISTRIBUTE
36 IF (lanafile) THEN
37#else
38 IF (lanafile.and.(tile.eq.0)) THEN
39#endif
40 ananame(11)=myfile
41 END IF
42!
43 RETURN
44 END SUBROUTINE ana_m2clima
45!
46!***********************************************************************
47 SUBROUTINE ana_m2clima_tile (ng, tile, model, &
48 & LBi, UBi, LBj, UBj, &
49 & IminS, ImaxS, JminS, JmaxS)
50!***********************************************************************
51!
52 USE mod_param
53 USE mod_clima
54 USE mod_scalars
55!
57#ifdef DISTRIBUTE
59#endif
60!
61! Imported variable declarations.
62!
63 integer, intent(in) :: ng, tile, model
64 integer, intent(in) :: LBi, UBi, LBj, UBj
65 integer, intent(in) :: IminS, ImaxS, JminS, JmaxS
66!
67! Local variable declarations.
68!
69 integer :: i, j
70
71#include "set_bounds.h"
72!
73!-----------------------------------------------------------------------
74! Set 2D momentum climatology.
75!-----------------------------------------------------------------------
76!
77 IF (lm2clm(ng)) THEN
78 DO j=jstrt,jendt
79 DO i=istrp,iendt
80 clima(ng)%ubarclm(i,j)=???
81 END DO
82 END DO
83 DO j=jstrp,jendt
84 DO i=istrt,iendt
85 clima(ng)%vbarclm(i,j)=???
86 END DO
87 END DO
88!
89! Exchange boundary data.
90!
91 IF (ewperiodic(ng).or.nsperiodic(ng)) THEN
92 CALL exchange_u2d_tile (ng, tile, &
93 & lbi, ubi, lbj, ubj, &
94 & clima(ng) % ubarclm)
95 CALL exchange_v2d_tile (ng, tile, &
96 & lbi, ubi, lbj, ubj, &
97 & clima(ng) % vbarclm)
98 END IF
99
100#ifdef DISTRIBUTE
101 CALL mp_exchange2d (ng, tile, model, 2, &
102 & lbi, ubi, lbj, ubj, &
103 & nghostpoints, &
104 & ewperiodic(ng), nsperiodic(ng), &
105 & clima(ng) % ubarclm, &
106 % CLIMA(ng) % vbarclm)
107#endif
108 END IF
109!
110 RETURN
111 END SUBROUTINE ana_m2clima_tile
subroutine ana_m2clima_tile(ng, tile, model, lbi, ubi, lbj, ubj, imins, imaxs, jmins, jmaxs)
Definition ana_m2clima.h:50
subroutine ana_m2clima(ng, tile, model)
Definition ana_m2clima.h:3
subroutine exchange_u2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
subroutine exchange_v2d_tile(ng, tile, lbi, ubi, lbj, ubj, a)
type(t_clima), dimension(:), allocatable clima
Definition mod_clima.F:153
logical lanafile
character(len=256), dimension(39) ananame
integer nghostpoints
Definition mod_param.F:710
logical, dimension(:), allocatable ewperiodic
logical, dimension(:), allocatable nsperiodic
logical, dimension(:), allocatable lm2clm
subroutine mp_exchange2d(ng, tile, model, nvar, lbi, ubi, lbj, ubj, nghost, ew_periodic, ns_periodic, a, b, c, d)