/[MITgcm]/MITgcm/pkg/exf/exf_mapfields.F
ViewVC logotype

Annotation of /MITgcm/pkg/exf/exf_mapfields.F

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


Revision 1.5 - (hide annotations) (download)
Tue Feb 18 05:33:54 2003 UTC (21 years, 3 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint50c_post, checkpoint50c_pre, checkpoint48i_post, checkpoint50, checkpoint50d_post, checkpoint50b_pre, checkpoint48f_post, checkpoint48h_post, checkpoint50a_post, checkpoint50f_pre, checkpoint50e_pre, checkpoint50e_post, checkpoint50d_pre, checkpoint49, checkpoint48g_post, checkpoint50b_post
Changes since 1.4: +24 -41 lines
Merging from release1_p12:
o Modifications for using pkg/exf with pkg/seaice
  - improved description of the various forcing configurations
  - added basic radiation bulk formulae to pkg/exf
  - units/sign fix for evap computation in exf_getffields.F
  - updated verification/global_with_exf/results/output.txt
o Added pkg/sbo for computing IERS Special Bureau for the Oceans
  (SBO) core products, including oceanic mass, center-of-mass,
  angular, and bottom pressure (see pkg/sbo/README.sbo).
o Lower bound for viscosity/diffusivity in pkg/kpp/kpp_routines.F
  to avoid negative values in shallow regions.
  - updated verification/natl_box/results/output.txt
  - updated verification/lab_sea/results/output.txt
o MPI gather, scatter: eesupp/src/gather_2d.F and scatter_2d.F
o Added useSingleCpuIO option (see PARAMS.h).
o Updated useSingleCpuIO option in mdsio_writefield.F to
  work with multi-field files, e.g., for single-file pickup.
o pkg/seaice:
  - bug fix in growth.F: QNET for no shortwave case
  - added HeffFile for specifying initial sea-ice thickness
  - changed SEAICE_EXTERNAL_FLUXES wind stress implementation
o Added missing /* */ to CPP comments in pkg/seaice, pkg/exf,
  kpp_transport_t.F, forward_step.F, and the_main_loop.F
o pkg/seaice:
  - adjoint-friendly modifications
  - added a SEAICE_WRITE_PICKUP at end of the_model_main.F

1 dimitri 1.5 c $Header: /u/gcmpack/MITgcm/pkg/exf/exf_mapfields.F,v 1.4 2002/12/28 10:11:11 dimitri Exp $
2 heimbach 1.1
3     #include "EXF_CPPOPTIONS.h"
4    
5    
6 heimbach 1.2 subroutine exf_mapfields( mythid )
7 heimbach 1.1
8     c ==================================================================
9 heimbach 1.2 c SUBROUTINE exf_mapfields
10 heimbach 1.1 c ==================================================================
11     c
12 dimitri 1.5 c o Map external forcing fields (ustress, vstress, hflux, sflux,
13     c swflux, apressure, climsss, climsst, etc.) onto ocean model
14     c arrays (fu, fv, Qnet, EmPmR, Qsw, pload, sss, sst, etc.).
15     c This routine is included to separate the ocean state estimation
16     c tool as much as possible from the ocean model. Unit and sign
17     c conventions can be customized using variables exf_outscal_*,
18     c which are set in exf_readparms.F. See the header files
19     c exf_fields.h and FFIELDS.h for definitions of the various input
20     c and output fields and for default unit and sign convetions.
21 heimbach 1.1 c
22     c started: Christian Eckert eckert@mit.edu 09-Aug-1999
23     c
24     c changed: Christian Eckert eckert@mit.edu 11-Jan-2000
25     c - Restructured the code in order to create a package
26     c for the MITgcmUV.
27     c
28     c Christian Eckert eckert@mit.edu 12-Feb-2000
29     c - Changed Routine names (package prefix: exf_)
30     c
31     c Patrick Heimbach, heimbach@mit.edu 06-May-2000
32     c - added and changed CPP flag structure for
33     c ALLOW_BULKFORMULAE, ALLOW_ATM_TEMP
34     c
35     c Patrick Heimbach, heimbach@mit.edu 23-May-2000
36     c - sign change of ustress/vstress incorporated into
37 dimitri 1.5 c scaling factors exf_outscal_ust, exf_outscal_vst
38 dimitri 1.4 c
39 dimitri 1.5 c mods for pkg/seaice: menemenlis@jpl.nasa.gov 20-Dec-2002
40 dimitri 1.4 c
41 heimbach 1.1 c ==================================================================
42 heimbach 1.2 c SUBROUTINE exf_mapfields
43 heimbach 1.1 c ==================================================================
44    
45     implicit none
46    
47     c == global variables ==
48    
49     #include "EEPARAMS.h"
50     #include "SIZE.h"
51     #include "FFIELDS.h"
52 heimbach 1.2 #include "exf_param.h"
53 heimbach 1.1 #include "exf_constants.h"
54     #include "exf_fields.h"
55     #include "exf_clim_fields.h"
56 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
57     # include "tamc.h"
58     # include "tamc_keys.h"
59     #endif
60 heimbach 1.1 c == routine arguments ==
61    
62     c mythid - thread number for this instance of the routine.
63    
64     integer mythid
65    
66     c == local variables ==
67    
68     integer bi,bj
69     integer i,j
70     integer jtlo
71     integer jthi
72     integer itlo
73     integer ithi
74     integer jmin
75     integer jmax
76     integer imin
77     integer imax
78    
79     c == end of interface ==
80    
81     jtlo = mybylo(mythid)
82     jthi = mybyhi(mythid)
83     itlo = mybxlo(mythid)
84     ithi = mybxhi(mythid)
85     jmin = 1-oly
86     jmax = sny+oly
87     imin = 1-olx
88     imax = snx+olx
89    
90     do bj = jtlo,jthi
91     do bi = itlo,ithi
92 heimbach 1.2
93     #ifdef ALLOW_AUTODIFF_TAMC
94     act1 = bi - myBxLo(myThid)
95     max1 = myBxHi(myThid) - myBxLo(myThid) + 1
96     act2 = bj - myByLo(myThid)
97     max2 = myByHi(myThid) - myByLo(myThid) + 1
98     act3 = myThid - 1
99     max3 = nTx*nTy
100     act4 = ikey_dynamics - 1
101     ikey = (act1 + 1) + act2*max1
102     & + act3*max1*max2
103     & + act4*max1*max2*max3
104     #endif /* ALLOW_AUTODIFF_TAMC */
105    
106 heimbach 1.1 do j = jmin,jmax
107     do i = imin,imax
108 dimitri 1.5 c Heat flux.
109     qnet(i,j,bi,bj) = exf_outscal_hfl*hflux(i,j,bi,bj)
110 heimbach 1.2 enddo
111     enddo
112 heimbach 1.1
113 heimbach 1.2
114     do j = jmin,jmax
115     do i = imin,imax
116 dimitri 1.5 c Salt flux.
117     empmr(i,j,bi,bj)= exf_outscal_sfl*sflux(i,j,bi,bj)
118 heimbach 1.2 enddo
119     enddo
120 heimbach 1.1
121 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
122     CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
123     #endif
124     do j = jmin,jmax
125     do i = imin,imax
126 heimbach 1.1 c Zonal wind stress.
127 heimbach 1.2 if (ustress(i,j,bi,bj).gt.2.0D0) then
128     ustress(i,j,bi,bj)=2.0D0
129     endif
130     enddo
131     enddo
132     #ifdef ALLOW_AUTODIFF_TAMC
133     CADJ STORE ustress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
134     #endif
135     do j = jmin,jmax
136     do i = imin,imax
137     if (ustress(i,j,bi,bj).lt.-2.0D0) then
138     ustress(i,j,bi,bj)=-2.0D0
139     endif
140     enddo
141     enddo
142     do j = jmin,jmax
143     do i = imin,imax
144 dimitri 1.5 fu(i,j,bi,bj) = exf_outscal_ust*ustress(i,j,bi,bj)
145 heimbach 1.2 enddo
146     enddo
147 heimbach 1.1
148 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
149     CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
150     #endif
151     do j = jmin,jmax
152     do i = imin,imax
153 heimbach 1.1 c Meridional wind stress.
154 heimbach 1.2 if (vstress(i,j,bi,bj).gt.2.0D0) then
155     vstress(i,j,bi,bj)=2.0D0
156     endif
157     enddo
158     enddo
159     #ifdef ALLOW_AUTODIFF_TAMC
160     CADJ STORE vstress(:,:,bi,bj) = comlev1_bibj, key=ikey, byte=isbyte
161     #endif
162     do j = jmin,jmax
163     do i = imin,imax
164     if (vstress(i,j,bi,bj).lt.-2.0D0) then
165     vstress(i,j,bi,bj)=-2.0D0
166     endif
167     enddo
168     enddo
169     do j = jmin,jmax
170     do i = imin,imax
171 dimitri 1.5 fv(i,j,bi,bj) = exf_outscal_vst*vstress(i,j,bi,bj)
172 heimbach 1.2 enddo
173     enddo
174 heimbach 1.1
175 dimitri 1.5 #ifdef SHORTWAVE_HEATING
176 heimbach 1.1 c Short wave radiative flux.
177 heimbach 1.2 do j = jmin,jmax
178     do i = imin,imax
179 dimitri 1.5 qsw(i,j,bi,bj) = exf_outscal_swf*swflux(i,j,bi,bj)
180 heimbach 1.2 enddo
181     enddo
182 heimbach 1.1 #endif
183    
184     #ifdef ALLOW_CLIMSST_RELAXATION
185 heimbach 1.2 do j = jmin,jmax
186     do i = imin,imax
187 dimitri 1.5 sst(i,j,bi,bj) = exf_outscal_sst*climsst(i,j,bi,bj)
188 heimbach 1.2 enddo
189     enddo
190 heimbach 1.1 #endif
191    
192     #ifdef ALLOW_CLIMSSS_RELAXATION
193 heimbach 1.2 do j = jmin,jmax
194     do i = imin,imax
195 dimitri 1.5 sss(i,j,bi,bj) = exf_outscal_sss*climsss(i,j,bi,bj)
196 heimbach 1.2 enddo
197     enddo
198 heimbach 1.1 #endif
199    
200 heimbach 1.2 #ifdef ATMOSPHERIC_LOADING
201     do j = jmin,jmax
202     do i = imin,imax
203 dimitri 1.5 pload(i,j,bi,bj)=exf_outscal_apressure*apressure(i,j,bi,bj)
204 heimbach 1.1 enddo
205     enddo
206 heimbach 1.2 #endif
207    
208 heimbach 1.1 enddo
209     enddo
210    
211     c Update the tile edges.
212    
213     _EXCH_XY_R4( qnet, mythid )
214     _EXCH_XY_R4( empmr, mythid )
215 cheisey 1.3 c _EXCH_XY_R4( fu, mythid )
216     c _EXCH_XY_R4( fv, mythid )
217     CALL EXCH_UV_XY_RS(fu, fv, .TRUE., myThid)
218 dimitri 1.5 #ifdef SHORTWAVE_HEATING
219 heimbach 1.1 _EXCH_XY_R4( qsw, mythid )
220     #endif
221     #ifdef ALLOW_CLIMSST_RELAXATION
222     _EXCH_XY_R4( sst, mythid )
223     #endif
224     #ifdef ALLOW_CLIMSSS_RELAXATION
225     _EXCH_XY_R4( sss, mythid )
226 heimbach 1.2 #endif
227     #ifdef ATMOSPHERIC_LOADING
228     _EXCH_XY_R4( pload, mythid )
229 heimbach 1.1 #endif
230    
231     end

  ViewVC Help
Powered by ViewVC 1.1.22