Changeset 3991 for branches/2013/dev_r3987_METO1_MERCATOR6_OBC_BDY_merge/NEMOGCM/NEMO/OPA_SRC/BDY/bdy_oce.F90
- Timestamp:
- 2013-07-29T11:04:44+02:00 (11 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
branches/2013/dev_r3987_METO1_MERCATOR6_OBC_BDY_merge/NEMOGCM/NEMO/OPA_SRC/BDY/bdy_oce.F90
r3651 r3991 29 29 REAL , POINTER, DIMENSION(:,:) :: nbw 30 30 REAL , POINTER, DIMENSION(:,:) :: nbd 31 REAL , POINTER, DIMENSION(:) :: flagu 32 REAL , POINTER, DIMENSION(:) :: flagv 31 REAL , POINTER, DIMENSION(:,:) :: nbdout 32 REAL , POINTER, DIMENSION(:,:) :: flagu 33 REAL , POINTER, DIMENSION(:,:) :: flagv 33 34 END TYPE OBC_INDEX 34 35 36 !! Logicals in OBC_DATA structure are true if the chosen algorithm requires this 37 !! field as external data. If true the data can come from external files 38 !! or model initial conditions. If false then no "external" data array 39 !! is required for this field. 40 35 41 TYPE, PUBLIC :: OBC_DATA !: Storage for external data 42 INTEGER, DIMENSION(2) :: nread 43 LOGICAL :: ll_ssh 44 LOGICAL :: ll_u2d 45 LOGICAL :: ll_v2d 46 LOGICAL :: ll_u3d 47 LOGICAL :: ll_v3d 48 LOGICAL :: ll_tem 49 LOGICAL :: ll_sal 36 50 REAL, POINTER, DIMENSION(:) :: ssh 37 51 REAL, POINTER, DIMENSION(:) :: u2d … … 42 56 REAL, POINTER, DIMENSION(:,:) :: sal 43 57 #if defined key_lim2 58 LOGICAL :: ll_frld 59 LOGICAL :: ll_hicif 60 LOGICAL :: ll_hsnif 44 61 REAL, POINTER, DIMENSION(:) :: frld 45 62 REAL, POINTER, DIMENSION(:) :: hicif … … 63 80 INTEGER :: nn_volctl !: = 0 the total volume will have the variability of the surface Flux E-P 64 81 ! ! = 1 the volume will be constant during all the integration. 65 INTEGER, DIMENSION(jp_bdy) :: nn_dyn2d! Choice of boundary condition for barotropic variables (U,V,SSH)66 INTEGER, DIMENSION(jp_bdy) :: nn_dyn2d_dta!: = 0 use the initial state as bdy dta ;82 CHARACTER(len=20), DIMENSION(jp_bdy) :: cn_dyn2d ! Choice of boundary condition for barotropic variables (U,V,SSH) 83 INTEGER, DIMENSION(jp_bdy) :: nn_dyn2d_dta !: = 0 use the initial state as bdy dta ; 67 84 !: = 1 read it in a NetCDF file 68 85 !: = 2 read tidal harmonic forcing from a NetCDF file 69 86 !: = 3 read external data AND tidal harmonic forcing from NetCDF files 70 INTEGER, DIMENSION(jp_bdy) :: nn_dyn3d! Choice of boundary condition for baroclinic velocities71 INTEGER, DIMENSION(jp_bdy) :: nn_dyn3d_dta!: = 0 use the initial state as bdy dta ;87 CHARACTER(len=20), DIMENSION(jp_bdy) :: cn_dyn3d ! Choice of boundary condition for baroclinic velocities 88 INTEGER, DIMENSION(jp_bdy) :: nn_dyn3d_dta !: = 0 use the initial state as bdy dta ; 72 89 !: = 1 read it in a NetCDF file 73 INTEGER, DIMENSION(jp_bdy) :: nn_tra! Choice of boundary condition for active tracers (T and S)74 INTEGER, DIMENSION(jp_bdy) :: nn_tra_dta!: = 0 use the initial state as bdy dta ;90 CHARACTER(len=20), DIMENSION(jp_bdy) :: cn_tra ! Choice of boundary condition for active tracers (T and S) 91 INTEGER, DIMENSION(jp_bdy) :: nn_tra_dta !: = 0 use the initial state as bdy dta ; 75 92 !: = 1 read it in a NetCDF file 76 93 LOGICAL, DIMENSION(jp_bdy) :: ln_tra_dmp !: =T Tracer damping 77 94 LOGICAL, DIMENSION(jp_bdy) :: ln_dyn3d_dmp !: =T Baroclinic velocity damping 78 95 REAL, DIMENSION(jp_bdy) :: rn_time_dmp !: Damping time scale in days 96 REAL, DIMENSION(jp_bdy) :: rn_time_dmp_out !: Damping time scale in days at radiation outflow points 79 97 80 98 #if defined key_lim2 81 INTEGER, DIMENSION(jp_bdy) :: nn_ice_lim2! Choice of boundary condition for sea ice variables82 INTEGER, DIMENSION(jp_bdy) :: nn_ice_lim2_dta!: = 0 use the initial state as bdy dta ;83 !: = 1 read it in a NetCDF file99 CHARACTER(len=20), DIMENSION(jp_bdy) :: nn_ice_lim2 ! Choice of boundary condition for sea ice variables 100 INTEGER, DIMENSION(jp_bdy) :: nn_ice_lim2_dta !: = 0 use the initial state as bdy dta ; 101 !: = 1 read it in a NetCDF file 84 102 #endif 85 103 ! … … 88 106 !! Global variables 89 107 !!---------------------------------------------------------------------- 90 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: bdytmask !: Mask defining computational domain at T-points 91 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: bdyumask !: Mask defining computational domain at U-points 92 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: bdyvmask !: Mask defining computational domain at V-points 108 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:), TARGET :: bdytmask !: Mask defining computational domain at T-points 109 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:), TARGET :: bdyumask !: Mask defining computational domain at U-points 110 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:), TARGET :: bdyvmask !: Mask defining computational domain at V-points 111 REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:), TARGET :: bdyfmask !: Mask defining computational domain at F-points 93 112 94 113 REAL(wp) :: bdysurftot !: Lateral surface of unstructured open boundary 95 114 96 REAL(wp), POINTER, DIMENSION(:,:) :: pssh !:97 REAL(wp), POINTER, DIMENSION(:,:) :: phur !:98 REAL(wp), POINTER, DIMENSION(:,:) :: phvr !: Pointers for barotropic fields99 REAL(wp), POINTER, DIMENSION(:,:) :: pu 2d!:100 REAL(wp), POINTER, DIMENSION(:,:) :: pv 2d!:115 REAL(wp), POINTER, DIMENSION(:,:) :: pssh !: 116 REAL(wp), POINTER, DIMENSION(:,:) :: phur !: 117 REAL(wp), POINTER, DIMENSION(:,:) :: phvr !: Pointers for barotropic fields 118 REAL(wp), POINTER, DIMENSION(:,:) :: pub2d, pun2d, pua2d !: 119 REAL(wp), POINTER, DIMENSION(:,:) :: pvb2d, pvn2d, pva2d !: 101 120 102 121 !!---------------------------------------------------------------------- … … 109 128 REAL(wp), ALLOCATABLE, DIMENSION(:,:,:), TARGET :: dta_global2 !: workspace for reading in global data arrays (struct. bdy) 110 129 TYPE(OBC_INDEX), DIMENSION(jp_bdy), TARGET :: idx_bdy !: bdy indices (local process) 111 TYPE(OBC_DATA) , DIMENSION(jp_bdy) 130 TYPE(OBC_DATA) , DIMENSION(jp_bdy), TARGET :: dta_bdy !: bdy external data (local process) 112 131 113 132 !!---------------------------------------------------------------------- … … 125 144 !!---------------------------------------------------------------------- 126 145 ! 127 ALLOCATE( bdytmask(jpi,jpj) , bdyumask(jpi,jpj), bdyvmask(jpi,jpj), 146 ALLOCATE( bdytmask(jpi,jpj) , bdyumask(jpi,jpj), bdyvmask(jpi,jpj), bdyfmask(jpi,jpj), & 128 147 & STAT=bdy_oce_alloc ) 129 148 !
Note: See TracChangeset
for help on using the changeset viewer.