/[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.1 - (hide annotations) (download)
Mon May 14 22:08:41 2001 UTC (23 years ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint44e_post, checkpoint46l_post, checkpoint46g_pre, checkpoint46f_post, checkpoint44f_post, checkpoint46b_post, checkpoint43a-release1mods, checkpoint40pre3, checkpoint40pre1, checkpoint40pre7, checkpoint40pre6, checkpoint40pre9, checkpoint40pre8, checkpoint46l_pre, chkpt44d_post, release1_p1, release1_p2, release1_p3, release1_p4, release1_p5, release1_p6, checkpoint44e_pre, release1_b1, checkpoint43, release1_chkpt44d_post, checkpoint46d_pre, checkpoint40pre2, release1-branch_tutorials, checkpoint45d_post, checkpoint46j_pre, chkpt44a_post, checkpoint44h_pre, checkpoint40pre4, checkpoint46a_post, checkpoint46j_post, checkpoint46k_post, chkpt44c_pre, checkpoint45a_post, ecco_c44_e19, ecco_c44_e18, ecco_c44_e17, ecco_c44_e16, checkpoint44g_post, checkpoint46e_pre, checkpoint45b_post, checkpoint46b_pre, release1-branch-end, release1_final_v1, checkpoint46c_pre, checkpoint46, checkpoint44b_post, checkpoint46h_pre, checkpoint46m_post, checkpoint46a_pre, checkpoint45c_post, checkpoint44h_post, checkpoint46g_post, checkpoint39, checkpoint40pre5, chkpt44a_pre, checkpoint46i_post, ecco_c44_e20, ecco_c44_e21, checkpoint46c_post, ecco-branch-mod1, ecco-branch-mod2, ecco-branch-mod3, ecco-branch-mod4, ecco-branch-mod5, checkpoint46e_post, release1_beta1, checkpoint44b_pre, checkpoint42, checkpoint40, checkpoint41, checkpoint44, checkpoint45, checkpoint46h_post, chkpt44c_post, checkpoint44f_pre, checkpoint46d_post, release1-branch_branchpoint
Branch point for: release1_final, release1-branch, release1, ecco-branch, release1_coupled
Added external forcing package.
Not presently supported by mitgcm, i.e. disabled by default.

1 heimbach 1.1 c $Header: /u/gcmpack/development/heimbach/ecco_env/pkg/exf/exf_mapfields.F,v 1.5 2001/02/02 19:43:47 heimbach Exp $
2    
3     #include "EXF_CPPOPTIONS.h"
4    
5    
6     subroutine exf_MapFields(
7     I mythid
8     & )
9    
10     c ==================================================================
11     c SUBROUTINE exf_MapFields
12     c ==================================================================
13     c
14     c o Map the external forcing fields on the ocean model arrays. This
15     c routine is included to separate the ocean state estimation tool
16     c as much as possible from the ocean model. Unit conversion factors
17     c are to be set by the user.
18     c
19     c The units have to be such that the individual forcing record has
20     c units equal to [quantity/s]. See the header file *FFIELDS.h* of
21     c the MITgcmuv.
22     c
23     c Required units such that no scaling has to be applied:
24     c
25     c heat flux: input file W/m^2
26     c salt flux: input file m/s
27     c zonal wind stress: input file N/m^2
28     c merid. wind stress: input file N/m^2
29     c
30     c To allow for such unit conversions this routine contains scaling
31     c factors scal_quantity.
32     c
33     c started: Christian Eckert eckert@mit.edu 09-Aug-1999
34     c
35     c changed: Christian Eckert eckert@mit.edu 11-Jan-2000
36     c
37     c - Restructured the code in order to create a package
38     c for the MITgcmUV.
39     c
40     c Christian Eckert eckert@mit.edu 12-Feb-2000
41     c
42     c - Changed Routine names (package prefix: exf_)
43     c
44     c Patrick Heimbach, heimbach@mit.edu 06-May-2000
45     c
46     c - added and changed CPP flag structure for
47     c ALLOW_BULKFORMULAE, ALLOW_ATM_TEMP
48     c
49     c Patrick Heimbach, heimbach@mit.edu 23-May-2000
50     c
51     c - sign change of ustress/vstress incorporated into
52     c scaling factors scal_ust, scal_vst
53     c
54     c ==================================================================
55     c SUBROUTINE exf_MapFields
56     c ==================================================================
57    
58     implicit none
59    
60     c == global variables ==
61    
62     #include "EEPARAMS.h"
63     #include "SIZE.h"
64     #include "FFIELDS.h"
65     #include "exf_constants.h"
66     #include "exf_fields.h"
67     #include "exf_clim_fields.h"
68    
69     c == routine arguments ==
70    
71     c mythid - thread number for this instance of the routine.
72    
73     integer mythid
74    
75     c == local variables ==
76    
77     integer bi,bj
78     integer i,j
79     integer jtlo
80     integer jthi
81     integer itlo
82     integer ithi
83     integer jmin
84     integer jmax
85     integer imin
86     integer imax
87     _RL scal_hfl
88     _RL scal_ust
89     _RL scal_vst
90     _RL scal_swf
91     _RL scal_sst
92     _RL scal_sss
93     #if (defined (ALLOW_BULKFORMULAE) && defined (ALLOW_ATM_TEMP))
94     _RL scal_prc
95     #else
96     _RL scal_sfl
97     #endif
98    
99     c == end of interface ==
100    
101     jtlo = mybylo(mythid)
102     jthi = mybyhi(mythid)
103     itlo = mybxlo(mythid)
104     ithi = mybxhi(mythid)
105     jmin = 1-oly
106     jmax = sny+oly
107     imin = 1-olx
108     imax = snx+olx
109    
110     scal_hfl = 1. _d 0
111     scal_ust = -1. _d 0
112     scal_vst = -1. _d 0
113     scal_swf = 1. _d 0
114     scal_sst = 1. _d 0
115     scal_sss = 1. _d 0
116     #if (defined (ALLOW_BULKFORMULAE) && defined (ALLOW_ATM_TEMP))
117     scal_prc = 1. _d 0
118     #else
119     scal_sfl = 1. _d 0
120     #endif
121    
122     do bj = jtlo,jthi
123     do bi = itlo,ithi
124     do j = jmin,jmax
125     do i = imin,imax
126    
127     c Heat flux.
128     qnet(i,j,bi,bj) = scal_hfl*hflux(i,j,bi,bj)
129    
130     c Salt flux.
131     #if (defined (ALLOW_BULKFORMULAE) && defined (ALLOW_ATM_TEMP))
132     empmr(i,j,bi,bj) = scal_prc*precip(i,j,bi,bj)
133     #else
134     empmr(i,j,bi,bj) = scal_sfl*sflux(i,j,bi,bj)
135     #endif
136    
137     c Zonal wind stress.
138     fu(i,j,bi,bj) = scal_ust*ustress(i,j,bi,bj)
139    
140     c Meridional wind stress.
141     fv(i,j,bi,bj) = scal_vst*vstress(i,j,bi,bj)
142    
143     #ifdef ALLOW_KPP || (defined (ALLOW_BULKFORMULAE) && defined (ALLOW_ATM_TEMP)))
144     c Short wave radiative flux.
145     qsw(i,j,bi,bj) = scal_swf*swflux(i,j,bi,bj)
146     #endif
147    
148     #ifdef ALLOW_CLIMSST_RELAXATION
149     sst(i,j,bi,bj) = scal_sst*climsst(i,j,bi,bj)
150     #endif
151    
152     #ifdef ALLOW_CLIMSSS_RELAXATION
153     sss(i,j,bi,bj) = scal_sss*climsss(i,j,bi,bj)
154     #endif
155    
156     enddo
157     enddo
158     enddo
159     enddo
160    
161     c Update the tile edges.
162    
163     _EXCH_XY_R4( qnet, mythid )
164     _EXCH_XY_R4( empmr, mythid )
165     _EXCH_XY_R4( fu, mythid )
166     _EXCH_XY_R4( fv, mythid )
167     #ifdef ALLOW_KPP
168     _EXCH_XY_R4( qsw, mythid )
169     #endif
170     #ifdef ALLOW_CLIMSST_RELAXATION
171     _EXCH_XY_R4( sst, mythid )
172     #endif
173     #ifdef ALLOW_CLIMSSS_RELAXATION
174     _EXCH_XY_R4( sss, mythid )
175     #endif
176    
177     end

  ViewVC Help
Powered by ViewVC 1.1.22