/[MITgcm]/MITgcm/pkg/fizhi/fizhi_utils.F
ViewVC logotype

Annotation of /MITgcm/pkg/fizhi/fizhi_utils.F

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


Revision 1.1 - (hide annotations) (download)
Wed Jun 9 21:35:31 2004 UTC (20 years ago) by molod
Branch: MAIN
Some utility routines called all over fizhi

1 molod 1.1 SUBROUTINE AVGGRD(IGRD,CHFR,B,MXCHPS,NCHP,FRACG,A,IM,JM)
2    
3     implicit none
4     real zero
5     parameter (zero = 0.)
6     integer i,im,jm,mxchps,nchp
7     integer igrd(mxchps)
8     real A(IM,JM), B(MXCHPS), CHFR(MXCHPS), FRACG(IM,JM)
9    
10     IF(NCHP.GT.0) THEN
11     DO I=1,IM*JM
12     A(I,1) = ZERO
13     ENDDO
14    
15     DO I=1,NCHP
16     A(IGRD(I),1) = A(IGRD(I),1) + B(I)*CHFR(I)
17     ENDDO
18    
19     DO I=1,IM*JM
20     IF(FRACG(I,1).GT.ZERO) THEN
21     A(I,1) = A(I,1)/FRACG(I,1)
22     ENDIF
23     ENDDO
24     ELSE
25     PRINT *, 'ERROR IN AVGGRD'
26     CALL my_EXIT(1)
27     ENDIF
28    
29     RETURN
30     END
31     SUBROUTINE BMPGRD(IGRD,CHFR,B,MXCHPS,NCHP,A,IM,JM)
32    
33     implicit none
34     integer i,im,jm,mxchps,nchp
35     integer igrd(mxchps)
36     real A(IM,JM), B(MXCHPS), CHFR(MXCHPS)
37    
38     IF(NCHP.GT.0) THEN
39     DO I=1,NCHP
40     A(IGRD(I),1) = A(IGRD(I),1) + B(I)*CHFR(I)
41     ENDDO
42     ELSE
43     PRINT *, 'ERROR IN BMPGRD'
44     CALL my_EXIT(1)
45     ENDIF
46    
47     RETURN
48     END
49     function minval (q,im)
50     implicit none
51     integer im, i
52     real q(im), minval
53     minval = 1.e15
54     do i=1,im
55     if( q(i).lt.minval ) minval = q(i)
56     enddo
57     return
58     end
59     FUNCTION ERRF (ARG)
60     C***********************************************************************
61     C FUNCTION ERRF
62     C PURPOSE
63     C COMPUTES ERROR FUNCTION OF ARGUMENT
64     C USAGE
65     C CALLED BY TRBFLX
66     C DESCRIPTION OF PARAMETERS
67     C ARG - INPUTED ARGUMENT
68     C REMARKS:
69     C USED TO COMPUTE FRACTIONAL CLOUD COVER AND LIQUID WATER CONTENT
70     C FROM TURBULENCE STATISTICS
71     C **********************************************************************
72    
73     PARAMETER ( AA1 = 0.254829592 )
74     PARAMETER ( AA2 = -0.284496736 )
75     PARAMETER ( AA3 = 1.421413741 )
76     PARAMETER ( AA4 = -1.453152027 )
77     PARAMETER ( AA5 = 1.061405429 )
78     PARAMETER ( PP = 0.3275911 )
79     PARAMETER ( X2 = AA2 / AA1 )
80     PARAMETER ( X3 = AA3 / AA2 )
81     PARAMETER ( X4 = AA5 / AA3 )
82     PARAMETER ( X5 = AA5 / AA4 )
83    
84     ERRF = 1.
85     AARG=ABS(ARG)
86    
87     IF ( AARG .LT. 4.0 ) THEN
88     TT = 1./(1.+PP*AARG)
89     ERRF = 1. -
90     1 (AA1*TT*(1.+X2*TT*(1.+X3*TT*(1.+X4*TT*(1.+X5*TT)))))
91     2 * EXP(-AARG*AARG)
92     ENDIF
93    
94     IF ( ARG .LT. 0.0 ) ERRF = -ERRF
95    
96     RETURN
97     END
98    
99     SUBROUTINE P2KAP ( PKE,PE,IRN,LM,PK,T1 )
100     C***********************************************************************
101     C PURPOSE
102     C CALCULATE PHILLIPS P**KAPPA
103     C ARGUEMENTS
104     C PKE .... SIGMA EDGE-LEVEL PRESSURE**KAPPA
105     C PE ..... SIGMA EDGE-LEVEL PRESSURE
106     C IRN .... HORIZONTAL DIMENSION
107     C LM ..... VERTICAL DIMENSION
108     C PK ..... SIGMA LEVEL PRESSURE**KAPPA
109     C T1 ..... TEMPORARY WORK ARRAY
110     C
111     C***********************************************************************
112     C* GODDARD LABORATORY FOR ATMOSPHERES
113     C***********************************************************************
114     C
115     PARAMETER ( ONE = 1.0 )
116    
117     real PKE(IRN,LM+1)
118     real PE(IRN,LM+1)
119     real PK(IRN,LM)
120     real T1(IRN,LM+1)
121    
122     AKAP = GETCON('KAPPA')
123    
124     IRNLM = IRN*LM
125     IRNLMP1 = IRN*(LM+1)
126    
127     DO I=1,IRNLMP1
128     T1(I,1) = PE(I,1)*PKE(I,1)
129     ENDDO
130    
131     DO I=1,IRNLM
132     PK(I,1) = ( T1(I,2)-T1(I,1) ) / ((ONE+AKAP)*( PE(I,2)-PE(I,1) ))
133     ENDDO
134    
135     RETURN
136     END
137     SUBROUTINE STRIP(A,B,IA,IB,L,K)
138     real A(IA,L), B(IB,L)
139     INTEGER OFFSET
140    
141     OFFSET = IB*(K-1)
142     LEN = MIN(IB,IA-OFFSET)
143     OFFSET = OFFSET+1
144    
145     IF(LEN.EQ.IB) THEN
146     DO 100 J=1,L
147     DO 100 I=1,LEN
148     B(I,J) = A(I+OFFSET-1,J)
149     100 CONTINUE
150     ELSE
151     DO 200 J=1,L
152     DO 300 I=1,LEN
153     B(I,J) = A(I+OFFSET-1,J)
154     300 CONTINUE
155     DO 400 I=1,IB-LEN
156     B(LEN+I,J) = A(LEN+OFFSET-1,J)
157     400 CONTINUE
158     200 CONTINUE
159     ENDIF
160    
161     RETURN
162     END
163     SUBROUTINE PASTE(B,A,IB,IA,L,K)
164     real A(IA,L), B(IB,L)
165     INTEGER OFFSET
166    
167     OFFSET = IB*(K-1)
168     LEN = MIN(IB,IA-OFFSET)
169     OFFSET = OFFSET+1
170    
171     DO 100 J=1,L
172     DO 100 I=1,LEN
173     A(I+OFFSET-1,J) = B(I,J)
174     100 CONTINUE
175    
176     RETURN
177     END
178     SUBROUTINE PSTBMP(B,A,IB,IA,L,K)
179     real A(IA,L), B(IB,L)
180     INTEGER OFFSET
181    
182     OFFSET = IB*(K-1)
183     LEN = MIN(IB,IA-OFFSET)
184     OFFSET = OFFSET+1
185    
186     DO 100 J=1,L
187     DO 100 I=1,LEN
188     A(I+OFFSET-1,J) = A(I+OFFSET-1,J) + B(I,J)
189     100 CONTINUE
190     C
191     RETURN
192     END
193     SUBROUTINE STRINT(A,B,IA,IB,L,K)
194     INTEGER A(IA,L), B(IB,L)
195     INTEGER OFFSET
196    
197     OFFSET = IB*(K-1)
198     LEN = MIN(IB,IA-OFFSET)
199     OFFSET = OFFSET+1
200    
201     IF(LEN.EQ.IB) THEN
202     DO 100 J=1,L
203     DO 100 I=1,LEN
204     B(I,J) = A(I+OFFSET-1,J)
205     100 CONTINUE
206     ELSE
207     DO 200 J=1,L
208     DO 300 I=1,LEN
209     B(I,J) = A(I+OFFSET-1,J)
210     300 CONTINUE
211     DO 400 I=1,IB-LEN
212     B(LEN+I,J) = A(LEN+OFFSET-1,J)
213     400 CONTINUE
214     200 CONTINUE
215     ENDIF
216    
217     RETURN
218     END
219     SUBROUTINE QSAT (TT,P,Q,DQDT,LDQDT)
220     C***********************************************************************
221     C
222     C PURPOSE:
223     C ========
224     C Compute Saturation Specific Humidity
225     C
226     C INPUT:
227     C ======
228     C TT ......... Temperature (Kelvin)
229     C P .......... Pressure (mb)
230     C LDQDT ...... Logical Flag to compute QSAT Derivative
231     C
232     C OUTPUT:
233     C =======
234     C Q .......... Saturation Specific Humidity
235     C DQDT ....... Saturation Specific Humidity Derivative wrt Temperature
236     C
237     C
238     C***********************************************************************
239     C* GODDARD LABORATORY FOR ATMOSPHERES *
240     C***********************************************************************
241    
242     IMPLICIT NONE
243     REAL TT, P, Q, DQDT
244     LOGICAL LDQDT
245     REAL AIRMW, H2OMW
246    
247     PARAMETER ( AIRMW = 28.97 )
248     PARAMETER ( H2OMW = 18.01 )
249    
250     REAL ESFAC, ERFAC
251     PARAMETER ( ESFAC = H2OMW/AIRMW )
252     PARAMETER ( ERFAC = (1.0-ESFAC)/ESFAC )
253    
254     real aw0, aw1, aw2, aw3, aw4, aw5, aw6
255     real bw0, bw1, bw2, bw3, bw4, bw5, bw6
256     real ai0, ai1, ai2, ai3, ai4, ai5, ai6
257     real bi0, bi1, bi2, bi3, bi4, bi5, bi6
258    
259     real d0, d1, d2, d3, d4, d5, d6
260     real e0, e1, e2, e3, e4, e5, e6
261     real f0, f1, f2, f3, f4, f5, f6
262     real g0, g1, g2, g3, g4, g5, g6
263    
264     c ********************************************************
265     c *** Polynomial Coefficients WRT Water (Lowe, 1977) ****
266     c *** (Valid +50 C to -50 C) ****
267     c ********************************************************
268    
269     parameter ( aw0 = 6.107799961e+00 * esfac )
270     parameter ( aw1 = 4.436518521e-01 * esfac )
271     parameter ( aw2 = 1.428945805e-02 * esfac )
272     parameter ( aw3 = 2.650648471e-04 * esfac )
273     parameter ( aw4 = 3.031240396e-06 * esfac )
274     parameter ( aw5 = 2.034080948e-08 * esfac )
275     parameter ( aw6 = 6.136820929e-11 * esfac )
276    
277     parameter ( bw0 = +4.438099984e-01 * esfac )
278     parameter ( bw1 = +2.857002636e-02 * esfac )
279     parameter ( bw2 = +7.938054040e-04 * esfac )
280     parameter ( bw3 = +1.215215065e-05 * esfac )
281     parameter ( bw4 = +1.036561403e-07 * esfac )
282     parameter ( bw5 = +3.532421810e-10 * esfac )
283     parameter ( bw6 = -7.090244804e-13 * esfac )
284    
285    
286     c ********************************************************
287     c *** Polynomial Coefficients WRT Ice (Lowe, 1977) ****
288     c *** (Valid +0 C to -50 C) ****
289     c ********************************************************
290    
291     parameter ( ai0 = +6.109177956e+00 * esfac )
292     parameter ( ai1 = +5.034698970e-01 * esfac )
293     parameter ( ai2 = +1.886013408e-02 * esfac )
294     parameter ( ai3 = +4.176223716e-04 * esfac )
295     parameter ( ai4 = +5.824720280e-06 * esfac )
296     parameter ( ai5 = +4.838803174e-08 * esfac )
297     parameter ( ai6 = +1.838826904e-10 * esfac )
298    
299     parameter ( bi0 = +5.030305237e-01 * esfac )
300     parameter ( bi1 = +3.773255020e-02 * esfac )
301     parameter ( bi2 = +1.267995369e-03 * esfac )
302     parameter ( bi3 = +2.477563108e-05 * esfac )
303     parameter ( bi4 = +3.005693132e-07 * esfac )
304     parameter ( bi5 = +2.158542548e-09 * esfac )
305     parameter ( bi6 = +7.131097725e-12 * esfac )
306    
307    
308     c ********************************************************
309     c *** Polynomial Coefficients WRT Ice ****
310     c *** Starr and Cox (1985) (Valid -40 C to -70 C) ****
311     c ********************************************************
312    
313    
314     parameter ( d0 = 0.535098336e+01 * esfac )
315     parameter ( d1 = 0.401390832e+00 * esfac )
316     parameter ( d2 = 0.129690326e-01 * esfac )
317     parameter ( d3 = 0.230325039e-03 * esfac )
318     parameter ( d4 = 0.236279781e-05 * esfac )
319     parameter ( d5 = 0.132243858e-07 * esfac )
320     parameter ( d6 = 0.314296723e-10 * esfac )
321    
322     parameter ( e0 = 0.469290530e+00 * esfac )
323     parameter ( e1 = 0.333092511e-01 * esfac )
324     parameter ( e2 = 0.102164528e-02 * esfac )
325     parameter ( e3 = 0.172979242e-04 * esfac )
326     parameter ( e4 = 0.170017544e-06 * esfac )
327     parameter ( e5 = 0.916466531e-09 * esfac )
328     parameter ( e6 = 0.210844486e-11 * esfac )
329    
330    
331     c ********************************************************
332     c *** Polynomial Coefficients WRT Ice ****
333     c *** Starr and Cox (1985) (Valid -65 C to -95 C) ****
334     c ********************************************************
335    
336     parameter ( f0 = 0.298152339e+01 * esfac )
337     parameter ( f1 = 0.191372282e+00 * esfac )
338     parameter ( f2 = 0.517609116e-02 * esfac )
339     parameter ( f3 = 0.754129933e-04 * esfac )
340     parameter ( f4 = 0.623439266e-06 * esfac )
341     parameter ( f5 = 0.276961083e-08 * esfac )
342     parameter ( f6 = 0.516000335e-11 * esfac )
343    
344     parameter ( g0 = 0.312654072e+00 * esfac )
345     parameter ( g1 = 0.195789002e-01 * esfac )
346     parameter ( g2 = 0.517837908e-03 * esfac )
347     parameter ( g3 = 0.739410547e-05 * esfac )
348     parameter ( g4 = 0.600331350e-07 * esfac )
349     parameter ( g5 = 0.262430726e-09 * esfac )
350     parameter ( g6 = 0.481960676e-12 * esfac )
351    
352     REAL TMAX, TICE
353     PARAMETER ( TMAX=323.15, TICE=273.16)
354    
355     REAL T, D, W, QX, DQX
356     T = MIN(TT,TMAX) - TICE
357     DQX = 0.
358     QX = 0.
359    
360     c Fitting for temperatures above 0 degrees centigrade
361     c ---------------------------------------------------
362     if(t.gt.0.) then
363     qx = aw0+T*(aw1+T*(aw2+T*(aw3+T*(aw4+T*(aw5+T*aw6)))))
364     if (ldqdt) then
365     dqx = bw0+T*(bw1+T*(bw2+T*(bw3+T*(bw4+T*(bw5+T*bw6)))))
366     endif
367     endif
368    
369     c Fitting for temperatures between 0 and -40
370     c ------------------------------------------
371     if( t.le.0. .and. t.gt.-40.0 ) then
372     w = (40.0 + t)/40.0
373     qx = w *(aw0+T*(aw1+T*(aw2+T*(aw3+T*(aw4+T*(aw5+T*aw6))))))
374     . + (1.-w)*(ai0+T*(ai1+T*(ai2+T*(ai3+T*(ai4+T*(ai5+T*ai6))))))
375     if (ldqdt) then
376     dqx = w *(bw0+T*(bw1+T*(bw2+T*(bw3+T*(bw4+T*(bw5+T*bw6))))))
377     . + (1.-w)*(bi0+T*(bi1+T*(bi2+T*(bi3+T*(bi4+T*(bi5+T*bi6))))))
378     endif
379     endif
380    
381     c Fitting for temperatures between -40 and -70
382     c --------------------------------------------
383     if( t.le.-40.0 .and. t.ge.-70.0 ) then
384     qx = d0+T*(d1+T*(d2+T*(d3+T*(d4+T*(d5+T*d6)))))
385     if (ldqdt) then
386     dqx = e0+T*(e1+T*(e2+T*(e3+T*(e4+T*(e5+T*e6)))))
387     endif
388     endif
389    
390     c Fitting for temperatures less than -70
391     c --------------------------------------
392     if(t.lt.-70.0) then
393     qx = f0+t*(f1+t*(f2+t*(f3+t*(f4+t*(f5+t*f6)))))
394     if (ldqdt) then
395     dqx = g0+t*(g1+t*(g2+t*(g3+t*(g4+t*(g5+t*g6)))))
396     endif
397     endif
398    
399     c Compute Saturation Specific Humidity
400     c ------------------------------------
401     D = (P-ERFAC*QX)
402     IF(D.LT.0.) THEN
403     Q = 1.0
404     IF (LDQDT) DQDT = 0.
405     ELSE
406     D = 1.0 / D
407     Q = MIN(QX * D,1.0)
408     IF (LDQDT) DQDT = (1.0 + ERFAC*Q) * D * DQX
409     ENDIF
410     RETURN
411     END
412     subroutine vqsat (tt,p,q,dqdt,ldqdt,n)
413     implicit none
414     integer i,n
415     logical ldqdt
416     real tt(n), p(n), q(n), dqdt(n)
417     #if CRAY
418     #if f77
419     cfpp$ expand (qsat)
420     #endif
421     #if f90
422     !DIR$ inline always qsat
423     #endif
424     #endif
425     do i=1,n
426     call qsat ( tt(i),p(i),q(i),dqdt(i),ldqdt )
427     enddo
428     return
429     end
430    
431     subroutine stripit(a,b,irun,ia,ib,l,k)
432     implicit none
433     integer ia,ib,irun,l,k
434     real a(ia,l), b(ib,l)
435     integer i,j,len,offset
436    
437     offset = ib*(k-1)
438     len = min(ib,irun-offset)
439     offset = offset+1
440    
441     if(len.eq.ib) then
442     do 100 j=1,l
443     do 100 i=1,len
444     b(i,j) = a(i+offset-1,j)
445     100 continue
446     else
447     do 200 j=1,l
448     do 300 i=1,len
449     b(i,j) = a(i+offset-1,j)
450     300 continue
451     do 400 i=1,ib-len
452     b(len+i,j) = a(len+offset-1,j)
453     400 continue
454     200 continue
455     endif
456     return
457     end
458    
459     subroutine stripitint(a,b,irun,ia,ib,l,k)
460     implicit none
461     integer ia,ib,irun,l,k,a(ia,l),b(ib,l)
462     integer i,j,len,offset
463    
464     offset = ib*(k-1)
465     len = min(ib,irun-offset)
466     offset = offset+1
467    
468     if(len.eq.ib) then
469     do 100 j=1,l
470     do 100 i=1,len
471     b(i,j) = a(i+offset-1,j)
472     100 continue
473     else
474     do 200 j=1,l
475     do 300 i=1,len
476     b(i,j) = a(i+offset-1,j)
477     300 continue
478     do 400 i=1,ib-len
479     b(len+i,j) = a(len+offset-1,j)
480     400 continue
481     200 continue
482     endif
483     return
484     end
485    
486     subroutine pastit(b,a,ib,ia,irun,L,k)
487     implicit none
488     real a(ia,l), b(ib,l)
489     integer ib,ia,L,k,irun,len,offset
490     integer i,j
491    
492     offset = ib*(k-1)
493     len = min(ib,irun-offset)
494     offset = offset+1
495    
496     do 100 j=1,L
497     do 100 i=1,len
498     a(i+offset-1,j) = b(i,j)
499     100 continue
500     return
501     end
502    
503     subroutine pstbitint(b,a,ib,ia,irun,l,k)
504     implicit none
505     real a(ia,l)
506     integer b(ib,l)
507     integer ib,ia,L,k,irun,len,offset
508     integer i,j
509    
510     offset = ib*(k-1)
511     len = min(ib,irun-offset)
512     offset = offset+1
513    
514     do 100 j=1,L
515     do 100 i=1,len
516     a(i+offset-1,j) = a(i+offset-1,j) + float(b(i,j))
517     100 continue
518     return
519     end
520    
521    
522     subroutine pstbmpit(b,a,ib,ia,irun,l,k)
523     implicit none
524     real a(ia,l), b(ib,l)
525     integer ib,ia,L,k,irun,len,offset
526     integer i,j
527    
528     offset = ib*(k-1)
529     len = min(ib,irun-offset)
530     offset = offset+1
531    
532     do 100 j=1,L
533     do 100 i=1,len
534     a(i+offset-1,j) = a(i+offset-1,j) + b(i,j)
535     100 continue
536     return
537     end
538    
539     subroutine strip2tile(a,index,b,irun,ia,ib,levs,npeice)
540     c-----------------------------------------------------------------------
541     c subroutine strip2tile - extract one processors worth of grid points
542     c from a grid space array to a stripped tile
543     c space array
544     c
545     c input:
546     c a - array to be stripped FROM [ia,levs]
547     c index - array of horizontal indeces of grid points to convert to
548     c tile space
549     c irun - number of points in array a that need to be stripped
550     c ia - inner of dimension of source array
551     c ib - inner dimension of target array AND the number of points
552     c in the target array to be filled
553     c levs - number of vertical levels AND outer dimension of arrays
554     c npeice - the current strip number to be filled
555     c output:
556     c b - array to be filled, ie, one processors field [ib,levs]
557     c-----------------------------------------------------------------------
558     implicit none
559     integer ia,ib,irun,levs,npeice
560     real a(ia,levs), b(ib,levs)
561     integer index(irun)
562     integer i,k,len,offset
563    
564     offset = ib*(npeice-1)
565     len = min(ib,irun-offset)
566     offset = offset+1
567    
568     if(len.eq.ib) then
569     do 100 k=1,levs
570     do 100 i=1,len
571     b(i,k) = a(index(i+offset-1),k)
572     100 continue
573     else
574     do 200 k=1,levs
575     do 300 i=1,len
576     b(i,k) = a(index(i+offset-1),k)
577     300 continue
578     do 400 i=1,ib-len
579     b(len+i,k) = a(index(len+offset-1),k)
580     400 continue
581     200 continue
582     endif
583     return
584     end
585    
586     subroutine paste2grd_old(b,index,chfr,ib,numpts,a,ia,levs,npeice)
587     c-----------------------------------------------------------------------
588     c subroutine paste2grd - paste one processors worth of grid points
589     c from a stripped tile array to a grid
590     c space array
591     c
592     c input:
593     c b - array to be pasted back into grid space [ib,levs]
594     c index - array of horizontal indeces of grid points to convert to
595     c tile space[numpts]
596     c chfr - fractional area covered by the tile [ib]
597     c ib - inner dimension of source array AND number of points in
598     c array a that need to be pasted
599     c numpts - total number of points which were stripped
600     c ia - inner of dimension of target array
601     c levs - number of vertical levels AND outer dimension of arrays
602     c npeice - the current strip number to be filled
603     c output:
604     c a - grid space array to be filled [ia,levs]
605     c
606     c IMPORTANT NOTE:
607     c
608     c This routine will result in roundoff differences if called from
609     c within a parallel region.
610     c-----------------------------------------------------------------------
611    
612     implicit none
613     integer ia,ib,levs,numpts,npeice
614     integer index(numpts)
615     real a(ia,levs), b(ib,levs), chfr(ib)
616    
617     integer i,L,offset,len
618    
619     offset = ib*(npeice-1)
620     len = min(ib,numpts-offset)
621     offset = offset+1
622    
623     do L = 1,levs
624     do i=1,len
625     a(index(i+offset-1),L) = a(index(i+offset-1),L) + b(i,L)*chfr(i)
626     enddo
627     enddo
628     return
629     end
630     subroutine paste2grd (b,index,chfr,ib,numpts,a,ia,levs,npeice,check)
631     c-----------------------------------------------------------------------
632     c subroutine paste2grd - paste one processors worth of grid points
633     c from a stripped tile array to a grid
634     c space array
635     c
636     c input:
637     c b - array to be pasted back into grid space [ib,levs]
638     c index - array of horizontal indeces of grid points to convert to
639     c tile space[numpts]
640     c chfr - fractional area covered by the tile [ib]
641     c ib - inner dimension of source array AND number of points in
642     c array a that need to be pasted
643     c numpts - total number of points which were stripped
644     c ia - inner of dimension of target array
645     c levs - number of vertical levels AND outer dimension of arrays
646     c npeice - the current strip number to be filled
647     c check - logical to check for undefined values
648     c output:
649     c a - grid space array to be filled [ia,levs]
650     c
651     c IMPORTANT NOTE:
652     c
653     c This routine will result in roundoff differences if called from
654     c within a parallel region.
655     c-----------------------------------------------------------------------
656    
657     implicit none
658     integer ia,ib,levs,numpts,npeice
659     integer index(numpts)
660     real a(ia,levs), b(ib,levs), chfr(ib)
661     logical check
662    
663     integer i,L,offset,len
664     real undef,getcon
665    
666     offset = ib*(npeice-1)
667     len = min(ib,numpts-offset)
668     offset = offset+1
669    
670     if( check ) then
671     undef = getcon('UNDEF')
672     do L= 1,levs
673     do i= 1,len
674     if( a(index(i+offset-1),L).eq.undef .or. b(i,L).eq.undef ) then
675     a(index(i+offset-1),L) = undef
676     else
677     a(index(i+offset-1),L) = a(index(i+offset-1),L) + b(i,L)*chfr(i)
678     endif
679     enddo
680     enddo
681     else
682     do L= 1,levs
683     do i= 1,len
684     a(index(i+offset-1),L) = a(index(i+offset-1),L) + b(i,L)*chfr(i)
685     enddo
686     enddo
687     endif
688    
689     return
690     end
691     SUBROUTINE GRD2MSC(A,IM,JM,IGRD,B,MXCHPS,NCHP)
692    
693     implicit none
694     integer im,jm,mxchps,nchp
695     integer igrd(mxchps)
696     real A(IM,JM), B(MXCHPS)
697    
698     integer i
699    
700     IF(NCHP.GE.0) THEN
701     DO I=1,NCHP
702     B(I) = A(IGRD(I),1)
703     ENDDO
704     ELSE
705     PRINT *, 'ERROR IN GRD2MSC'
706     ENDIF
707    
708     RETURN
709     END
710    
711     SUBROUTINE MSC2GRD(IGRD,CHFR,B,MXCHPS,NCHP,FRACG,A,IM,JM)
712    
713     implicit none
714     real zero,one
715     parameter ( one = 1.)
716     parameter (zero = 0.)
717     integer im,jm,mxchps,nchp
718     integer igrd(mxchps)
719     real A(IM,JM), B(MXCHPS), CHFR(MXCHPS), FRACG(IM,JM)
720    
721     real VT1(IM,JM)
722     integer i
723    
724     IF(NCHP.GE.0) THEN
725     DO I=1,IM*JM
726     VT1(I,1) = ZERO
727     ENDDO
728    
729     DO I=1,NCHP
730     VT1(IGRD(I),1) = VT1(IGRD(I),1) + B(I)*CHFR(I)
731     ENDDO
732    
733     DO I=1,IM*JM
734     A(I,1) = A(I,1)*(ONE-FRACG(I,1)) + VT1(I,1)
735     ENDDO
736     ELSE
737     PRINT *, 'ERROR IN MSC2GRD'
738     ENDIF
739    
740     RETURN
741     END
742    
743     subroutine chpprm(nymd,nhms,mxchps,nchp,chlt,ityp,alai,
744     1 agrn,zoch,z2ch,cdrc,cdsc,sqsc,ufac,rsl1,rsl2,rdcs)
745    
746     implicit none
747    
748     integer nymd,nhms,nchp,mxchps,ityp(mxchps)
749     real chlt(mxchps)
750     real alai(mxchps),agrn(mxchps)
751     real zoch(mxchps), z2ch(mxchps), cdrc(mxchps), cdsc(mxchps)
752     real sqsc(mxchps), ufac(mxchps), rsl1(mxchps), rsl2(mxchps)
753     real rdcs(mxchps)
754    
755     C*********************************************************************
756     C********************* SUBROUTINE CHPPRM ****************************
757     C********************** 14 JUNE 1991 ******************************
758     C*********************************************************************
759    
760     integer ntyps
761     parameter (ntyps=10)
762    
763     real pblzet
764     parameter (pblzet = 50.)
765     integer k1,k2,nymd1,nhms1,nymd2,nhms2,i
766     real getcon,vkrm,rootl,vroot,dum1,dum2,alphaf
767     real facm,facp
768     real scat,d
769    
770     real
771     & vgdd(12, ntyps), vgz0(12, ntyps),
772     & vgrd(12, ntyps), vgrt(12, ntyps),
773    
774     & vgrf11(ntyps), vgrf12(ntyps),
775     & vgtr11(ntyps), vgtr12(ntyps),
776     & vgroca(ntyps), vgrotd(ntyps),
777     & vgrdrs(ntyps), vgz2 (ntyps)
778    
779    
780     data vgz0 /
781     1 2.6530, 2.6530, 2.6530, 2.6530, 2.6530, 2.6530, 2.6530,
782     1 2.6530, 2.6530, 2.6530, 2.6530, 2.6530,
783     2 0.5200, 0.5200, 0.6660, 0.9100, 1.0310, 1.0440, 1.0420,
784     2 1.0370, 1.0360, 0.9170, 0.6660, 0.5200,
785     3 1.1120, 1.1030, 1.0880, 1.0820, 1.0760, 1.0680, 1.0730,
786     3 1.0790, 1.0820, 1.0880, 1.1030, 1.1120,
787     4 0.0777, 0.0778, 0.0778, 0.0779, 0.0778, 0.0771, 0.0759,
788     4 0.0766, 0.0778, 0.0779, 0.0778, 0.0778,
789     5 0.2450, 0.2450, 0.2270, 0.2000, 0.2000, 0.2000, 0.2000,
790     5 0.267, 0.292, 0.280, 0.258, 0.2450,
791     6 0.0752, 0.0752, 0.0752, 0.0752, 0.0752, 0.0757, 0.0777,
792     6 0.0778, 0.0774, 0.0752, 0.0752, 0.0752,
793     7 0.0112, 0.0112, 0.0112, 0.0112, 0.0112, 0.0112, 0.0112,
794     7 0.0112, 0.0112, 0.0112, 0.0112, 0.0112,
795     8 0.0112, 0.0112, 0.0112, 0.0112, 0.0112, 0.0112, 0.0112,
796     8 0.0112, 0.0112, 0.0112, 0.0112, 0.0112,
797     9 0.0112, 0.0112, 0.0112, 0.0112, 0.0112, 0.0112, 0.0112,
798     9 0.0112, 0.0112, 0.0112, 0.0112, 0.0112,
799     1 0.0112, 0.0112, 0.0112, 0.0112, 0.0112, 0.0112, 0.0112,
800     1 0.0112, 0.0112, 0.0112, 0.0112, 0.0112
801     & /
802    
803    
804     data vgrd /
805     1 285.87, 285.87, 285.87, 285.87, 285.87, 285.87, 285.87,
806     1 285.87, 285.87, 285.87, 285.87, 285.87,
807     2 211.32, 211.32, 218.78, 243.40, 294.87, 345.90, 355.18,
808     2 341.84, 307.22, 244.84, 218.78, 211.32,
809     3 565.41, 587.05, 623.46, 638.13, 652.86, 675.04, 660.24,
810     3 645.49, 638.13, 623.46, 587.05, 565.41,
811     4 24.43, 24.63, 24.80, 24.96, 25.72, 27.74, 30.06,
812     4 28.86, 25.90, 25.11, 24.80, 24.63,
813     5 103.60, 103.60, 102.35, 100.72, 100.72, 100.72, 100.72,
814     5 105.30, 107.94, 106.59, 104.49, 103.60,
815     6 22.86, 22.86, 22.86, 22.86, 22.86, 23.01, 24.36,
816     6 24.69, 24.04, 22.86, 22.86, 22.86,
817     7 23.76, 23.76, 23.76, 23.76, 23.76, 23.76, 23.76,
818     7 23.76, 23.76, 23.76, 23.76, 23.76,
819     8 23.76, 23.76, 23.76, 23.76, 23.76, 23.76, 23.76,
820     8 23.76, 23.76, 23.76, 23.76, 23.76,
821     9 23.76, 23.76, 23.76, 23.76, 23.76, 23.76, 23.76,
822     9 23.76, 23.76, 23.76, 23.76, 23.76,
823     1 23.76, 23.76, 23.76, 23.76, 23.76, 23.76, 23.76,
824     1 23.76, 23.76, 23.76, 23.76, 23.76
825     & /
826    
827     data vgrt /
828     1 19737.8, 19737.8, 19737.8, 19737.8, 19737.8, 19737.8, 19737.8,
829     1 19737.8, 19737.8, 19737.8, 19737.8, 19737.8,
830     2 5010.0, 5010.0, 5270.0, 6200.0, 8000.0, 9700.0, 9500.0,
831     2 8400.0, 6250.0, 5270.0, 5010.0, 5010.0,
832     3 9000.0, 9200.0, 9533.3, 9666.7, 9800.0, 9866.7, 9733.3,
833     3 9666.7, 9533.3, 9200.0, 9000.0, 9000.0,
834     4 5500.0, 5625.0, 5750.0, 5875.0, 6625.0, 8750.0, 9375.0,
835     4 6875.0, 6000.0, 5750.0, 5625.0, 5500.0,
836     5 6500.0, 6000.0, 5500.0, 5500.0, 5500.0, 5500.0, 5500.0,
837     5 7500.0, 8500.0, 7000.0, 6500.0, 6500.0,
838     6 10625.0, 10625.0, 10625.0, 10625.0, 10625.0, 11250.0, 18750.0,
839     6 17500.0, 10625.0, 10625.0, 10625.0, 10625.0,
840     7 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
841     7 1.0, 1.0, 1.0, 1.0, 1.0,
842     8 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
843     8 1.0, 1.0, 1.0, 1.0, 1.0,
844     9 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
845     9 1.0, 1.0, 1.0, 1.0, 1.0,
846     1 1.0, 1.0, 1.0, 1.0, 1.0, 1.0, 1.0,
847     1 1.0, 1.0, 1.0, 1.0, 1.0
848     & /
849    
850    
851     data vgdd /
852     1 27.37, 27.37, 27.37, 27.37, 27.37, 27.37, 27.37,
853     1 27.37, 27.37, 27.37, 27.37, 27.37,
854     2 13.66, 13.66, 14.62, 15.70, 16.33, 16.62, 16.66,
855     2 16.60, 16.41, 15.73, 14.62, 13.66,
856     3 13.76, 13.80, 13.86, 13.88, 13.90, 13.93, 13.91,
857     3 13.89, 13.88, 13.86, 13.80, 13.76,
858     4 0.218, 0.227, 0.233, 0.239, 0.260, 0.299, 0.325,
859     4 0.313, 0.265, 0.244, 0.233, 0.227,
860     5 2.813, 2.813, 2.662, 2.391, 2.391, 2.391, 2.391,
861     5 2.975, 3.138, 3.062, 2.907, 2.813,
862     6 0.10629, 0.10629, 0.10629, 0.10629, 0.10629, 0.12299, 0.21521,
863     6 0.22897, 0.19961, 0.10629, 0.10629, 0.10629,
864     7 0.0001, 0.0001, 0.0001, 0.0001, 0.0001, 0.0001, 0.0001,
865     7 0.0001, 0.0001, 0.0001, 0.0001, 0.0001,
866     8 0.0001, 0.0001, 0.0001, 0.0001, 0.0001, 0.0001, 0.0001,
867     8 0.0001, 0.0001, 0.0001, 0.0001, 0.0001,
868     9 0.0001, 0.0001, 0.0001, 0.0001, 0.0001, 0.0001, 0.0001,
869     9 0.0001, 0.0001, 0.0001, 0.0001, 0.0001,
870     1 0.0001, 0.0001, 0.0001, 0.0001, 0.0001, 0.0001, 0.0001,
871     1 0.0001, 0.0001, 0.0001, 0.0001, 0.0001
872     & /
873    
874    
875     data vgrf11 /0.10,0.10,0.07,0.105,0.10,0.10,.001,.001,.001,.001/
876    
877     data vgrf12 /0.16,0.16,0.16,0.360,0.16,0.16,.001,.001,.001,.001/
878    
879     data vgtr11 /0.05,0.05,0.05,0.070,0.05,0.05,.001,.001,.001,.001/
880    
881     data vgtr12 /.001,.001,.001, .220,.001,.001,.001,.001,.001,.001/
882    
883     data vgroca /
884     & 0.384E-6, 0.384E-6, 0.384E-6, 0.384E-6, 0.384E-6, 0.384E-6,
885     & .1E-6, .1E-6, .1E-6, .1E-6 /
886    
887     data vgrotd /1.00,1.00,0.50,0.50,0.50,0.20,0.10,0.10,0.10,0.10/
888    
889     data vgrdrs /
890     & 0.75E13, 0.75E13, 0.75E13, 0.40E13, 0.75E13, 0.75E13,
891     & 0.10E13, 0.10E13, 0.10E13, 0.10E13 /
892    
893     data vgz2 /35.0, 20.0, 17.0, 0.6, 5.0, 0.6, 0.1, 0.1, 0.1, 0.1/
894    
895     vkrm = GETCON('VON KARMAN')
896    
897     call time_bound ( nymd,nhms, nymd1,nhms1, nymd2,nhms2, k1,k2 )
898     call interp_time ( nymd,nhms, nymd1,nhms1, nymd2,nhms2, facm,facp )
899    
900     do i=1,nchp
901    
902     zoch(i) = vgz0(k2,ityp(i))*facp + vgz0(k1,ityp(i))*facm
903     rdcs(i) = vgrd(k2,ityp(i))*facp + vgrd(k1,ityp(i))*facm
904    
905     rootl = vgrt(k2,ityp(i))*facp + vgrt(k1,ityp(i))*facm
906    
907     vroot = rootl * vgroca(ityp (i))
908     dum1 = log (vroot / (1. - vroot))
909     dum2 = 1. / (8. * 3.14159 * rootl)
910     alphaf = dum2 * (vroot - 3. -2. * dum1)
911    
912     rsl1(i) = vgrdrs (ityp (i)) / (rootl * vgrotd (ityp (i)))
913     rsl2(i) = alphaf / vgrotd (ityp (i))
914    
915     scat = agrn(i) *(vgtr11(ityp(i)) + vgrf11(ityp(i)))
916     & + (1. - agrn(i))*(vgtr12(ityp(i)) + vgrf12(ityp(i)))
917     sqsc(i) = sqrt (1. - scat)
918    
919     d = vgdd(k2,ityp(i))*facp + vgdd(k1,ityp(i))*facm
920     ufac(i) = log( (vgz2(ityp(i)) - d) / zoch(i) )
921     * / log( pblzet / zoch(i) )
922    
923     z2ch(i) = vgz2(ityp (i))
924    
925     cdsc(i) = pblzet/zoch(i)+1.
926     cdrc(i) = vkrm/log(cdsc(i))
927     cdrc(i) = cdrc(i)*cdrc(i)
928     cdsc(i) = sqrt(cdsc(i))
929     cdsc(i) = cdrc(i)*cdsc(i)
930    
931     enddo
932    
933     return
934     end

  ViewVC Help
Powered by ViewVC 1.1.22