source: trunk/LMDZ.MARS/libf/phymars/nlte_tcool.F @ 706

Last change on this file since 706 was 695, checked in by acolaitis, 12 years ago

GCM: modifications to allow for nesting compilation. Transparent to GCM user. MESOSCALE: handling of functions and modules for nesting compilation. Full nested configuration with new physics is now fully compiling.

File size: 50.0 KB
Line 
1c***********************************************************************
2                                                           
3      subroutine NLTEdlvr09_TCOOL (ngridgcm,n_gcm, 
4     @     p_gcm, t_gcm, z_gcm,
5     @     co2vmr_gcm, n2vmr_gcm, covmr_gcm, o3pvmr_gcm, 
6     @     q15umco2_gcm )
7
8c       jul 2011 malv+fgg                             
9c***********************************************************************
10                                                           
11      implicit none                                 
12         
13      include "dimensions.h"
14      include "dimphys.h"
15      include 'nlte_paramdef.h'
16      include 'nlte_commons.h'
17      include "chimiedata.h"
18      include "conc.h"
19
20c Arguments
21      integer n_gcm,ngridgcm
22      real p_gcm(ngridgcm,n_gcm), t_gcm(ngridgcm,n_gcm)
23      real co2vmr_gcm(ngridgcm,n_gcm), n2vmr_gcm(ngridgcm,n_gcm)
24      real covmr_gcm(ngridgcm,n_gcm), o3pvmr_gcm(ngridgcm,n_gcm)
25      real q15umco2_gcm(ngridgcm,n_gcm)
26      real z_gcm(ngridgcm,n_gcm)
27                                                             
28c local variables and constants                 
29      integer   iz, i, j, k, l, ig,istyle
30     
31      real*8            q15umco2_nl(nl)                       
32      real*8            zld(nl), zgcmd(n_gcm)                     
33      real*8          auxdgcm(n_gcm)                       
34
35
36      real p_ig(n_gcm),z_ig(n_gcm)
37      real t_ig(n_gcm)
38      real co2_ig(n_gcm),n2_ig(n_gcm),co_ig(n_gcm),o3p_ig(n_gcm)
39      real mmean_ig(n_gcm),cpnew_ig(n_gcm)
40       
41
42c**********************************************************************
43
44      do ig=1,ngridgcm
45         do l=1,n_gcm
46            p_ig(l)=p_gcm(ig,l)
47            t_ig(l)=t_gcm(ig,l)
48            co2_ig(l)=co2vmr_gcm(ig,l)
49            n2_ig(l)=n2vmr_gcm(ig,l)
50            o3p_ig(l)=o3pvmr_gcm(ig,l)
51            co_ig(l)=covmr_gcm(ig,l)
52            z_ig(l)=z_gcm(ig,l)/1000.
53            mmean_ig(l)=mmean(ig,l)
54            cpnew_ig(l)=cpnew(ig,l)
55         enddo
56
57         call NLTEdlvr09_ZGRID (n_gcm, 
58     @        p_ig, t_ig, z_ig,
59     @        co2_ig,n2_ig,co_ig,
60     $        o3p_ig , mmean_ig, cpnew_ig)
61
62c     And sets zero to all Curtis Matrixes and Escape Transmissions
63         call leetvt
64         call zero3m (c110,cup110,cdw110, nl)
65         call zero2v (taugr110,vc110, nl)
66         if (itt_cza.eq.24) then
67            call mzescape ( ig,taustar11,tauinf110,tauii110,
68     @           1, 1,irw_mztf,imu )
69            istyle=2
70            call mzescape_normaliz ( taustar11, istyle )   
71         else
72            call mztud (ig, c110,cup110,cdw110, vc110,taugr110,           
73     @           1, 1, irw_mztf, imu, 0,0,0 )
74         endif
75         call mztvc (ig,vc210, 1, 2, irw_mztf, imu, 0,0,0 )
76         call mztvc (ig,vc310, 1, 3, irw_mztf, imu, 0,0,0 )
77         call mztvc (ig,vc410, 1, 4, irw_mztf, imu, 0,0,0 )
78
79         call mzescape_fb (ig)       
80         input_cza = 0 
81         call NLTEdlvr09_CZALU(ig)
82         
83         if (itt_cza.ne.24) then
84            call mzescape_fh (ig)       
85            input_cza = 1 
86            call NLTEdlvr09_CZALU(ig)
87         endif
88         
89c total cooling rate                                               
90c smoothing and
91c interpolation back to original Pgrid
92c
93         do i = 1, nl                                   
94            q15umco2_nl(i) = hr110(i) + hr210(i) + hr310(i) + hr410(i)
95     @           + hr121(i)
96         enddo             
97         
98         do i=1,nl                                     
99            zld(i) = - dble ( alog(pl(i)) )                     
100         enddo                                         
101         do i=1,n_gcm                                     
102            zgcmd(i) = - dble( alog(p_gcm(ig,i)) ) 
103         enddo                   
104         call zerov( auxdgcm, n_gcm )
105         call interdp_limits                           
106     @        (auxdgcm,zgcmd,n_gcm,jlowerboundary,jtopboundary,
107     @        q15umco2_nl,zld,nl,1,nl,1)       
108         call suaviza ( auxdgcm, n_gcm, 1, zgcmd )     
109         
110         do i=1,n_gcm                                     
111            q15umco2_gcm(ig,i) = sngl( auxdgcm(i) )                       
112         enddo
113         
114      enddo
115       
116       
117c end subroutine                                   
118      return
119      end
120
121
122
123c***********************************************************************
124
125      subroutine NLTEdlvr09_ZGRID (n_gcm, 
126     @     p_gcm, t_gcm, z_gcm,
127     @     co2vmr_gcm, n2vmr_gcm, covmr_gcm, o3pvmr_gcm ,mmean_gcm,
128     @     cpnew_gcm)
129
130c     jul 2011 malv+fgg    First version
131c***********************************************************************
132                                               
133      implicit none                                 
134       
135      include "dimensions.h"
136      include "dimphys.h"
137      include 'nlte_paramdef.h'
138      include 'nlte_commons.h'
139      include 'chimiedata.h'
140      include 'conc.h'
141     
142c     Arguments
143      integer n_gcm
144      real p_gcm(n_gcm), t_gcm(n_gcm)
145      real co2vmr_gcm(n_gcm), n2vmr_gcm(n_gcm)
146      real covmr_gcm(n_gcm), o3pvmr_gcm(n_gcm)
147      real z_gcm(n_gcm)
148      real mmean_gcm(n_gcm)
149      real cpnew_gcm(n_gcm)
150
151c     local variables                               
152      integer i, j  , iz                               
153!     real  distancia, meanm, gz, Hkm
154      real  zmin, zmax, deltazz, deltazzy
155      real  nt_gcm(n_gcm)
156      real  mmean_nlte(n_gcm),cpnew_nlte(n_gcm)
157             
158c functions                                     
159      external  hrkday_convert                       
160      real              hrkday_convert                         
161                                               
162c***********************************************************************
163
164
165! Define working grid for MZ1D model (NL, ZL, ZMIN)
166! y otro mas fino para M.Curtis (NZ, ZX, ZXMIN = ZMIN
167
168! Para ello hace falta una z de ref del GCM, que voy a suponer la inferior
169
170! Primero, construimos escala z_gcm
171
172!       z_gcm (1) = zmin_gcm             ! [km]
173
174        !write (*,*) ' iz, p, g, H, z =', 1, p_gcm(1), z_gcm(1)
175!       do iz = 2, n_gcm
176!       do iz=1,n_gcm
177!          z_gcm(iz)=zlay(iz)/1.e3
178
179!         meanm = ( co2vmr_gcm(iz)*44. + o3pvmr_gcm(iz)*16.
180!     @               + n2vmr_gcm(iz)*28. + covmr_gcm(iz)*28. )
181!         meanm = meanm / n_avog
182!         distancia = ( radio + z_gcm(iz-1) )*1.e5
183!         gz = gg * masa / ( distancia * distancia )
184!          Hkm = 0.5*( t_gcm(iz)+t_gcm(iz-1) ) / ( meanm * gz )
185!          Hkm = kboltzman * Hkm *1e-5                           ! [km]
186!          z_gcm(iz) = z_gcm(iz-1) - Hkm * log( p_gcm(iz)/p_gcm(iz-1) )
187
188          !write (*,*) iz, p_gcm(iz), gz, Hkm, z_gcm(iz)
189
190!        enddo
191! Segundo, definimos los límites del modelo, entre las 2 presiones clave
192
193        ! Bottom boundary for NLTE model : Pbottom=2e-2mb=1.974e-5 atm
194      jlowerboundary = 1
195      do while ( p_gcm(jlowerboundary) .gt. Pbottom_atm )
196         jlowerboundary = jlowerboundary + 1
197      enddo
198      zmin = z_gcm(jlowerboundary)
199!       write (*,*) ' jlowerboundary, Pmin, zmin =',
200!     @            jlowerboundary, p_gcm(jlowerboundary), zmin
201
202        ! Top boundary for NLTE model : Ptop=2e-7mb = 1.974e-5 atm
203      jtopboundary = jlowerboundary 
204      do while ( p_gcm(jtopboundary) .gt. Ptop_atm )
205         jtopboundary = jtopboundary + 1
206      enddo
207      zmax = z_gcm(jtopboundary)
208!       write (*,*) ' jtopboundary, Pmax, zmax =',
209!     @      jtopboundary, p_gcm(jtopboundary),zmax
210
211      deltaz = (zmax-zmin) / (nl-1)
212      do i=1,nl                                     
213         zl(i) = zmin + (i-1) * deltaz
214      enddo                                         
215!       write (*,*) ' ZL grid:  dz,zmin,zmax ', deltaz, zl(1),zl(nl)
216! Creamos el perfil interpolando
217      call intersp (    pl,zl,nl,      p_gcm,z_gcm,n_gcm, 2) ! [atm]
218      call intersp (     t,zl,nl,      t_gcm,z_gcm,n_gcm, 1)       
219      do i = 1, n_gcm
220         nt_gcm(i) = 7.339e+21 * p_gcm(i) / t_gcm(i) ! [cm-3]     
221      enddo
222      call intersp (    nt,zl,nl,     nt_gcm,z_gcm,n_gcm, 2)       
223      call intersp (co2vmr,zl,nl, co2vmr_gcm,z_gcm,n_gcm, 1)       
224      call intersp ( n2vmr,zl,nl,  n2vmr_gcm,z_gcm,n_gcm, 1)       
225      call intersp ( covmr,zl,nl,  covmr_gcm,z_gcm,n_gcm, 1)       
226      call intersp (o3pvmr,zl,nl, o3pvmr_gcm,z_gcm,n_gcm, 1) 
227      call intersp (mmean_nlte,zl,nl,mmean_gcm,z_gcm,n_gcm,1)
228      call intersp (cpnew_nlte,zl,nl,cpnew_gcm,z_gcm,n_gcm,1)
229       
230
231      do i = 1, nl
232
233         co2(i) = nt(i) * co2vmr(i)
234         n2(i) = nt(i) * n2vmr(i)
235         co(i) = nt(i) * covmr(i)
236         o3p(i) = nt(i) * o3pvmr(i)
237
238!               hrkday_factor(i) =  hrkday_convert( t(i),       
239!     @           co2vmr(i), o3pvmr(i), n2vmr(i), covmr(i) )
240         hrkday_factor(i) = hrkday_convert(mmean_nlte(i)
241     &                          ,cpnew_nlte(i))
242
243      enddo
244                                               
245                                             
246
247c  Fine grid for transmittance calculations
248
249      deltazy = (zmax-zmin) / (nzy-1)
250      do i=1,nzy                                     
251         zy(i) = zmin + (i-1) * deltazy             
252      enddo                                         
253!       write (*,*) ' ZY grid:  nzy,dzy,zmin,zmax ',
254!     @         nzy, deltazy, zy(1),zy(nzy)
255
256      call intersp (    py,zy,nzy,      p_gcm,z_gcm,n_gcm, 2) ! [atm]
257      call intersp (    ty,zy,nzy,      t_gcm,z_gcm,n_gcm, 1)       
258      call intersp (   nty,zy,nzy,     nt_gcm,z_gcm,n_gcm, 2)       
259     
260      call intersp (  co2y,zy,nzy,   co2vmr_gcm,z_gcm,n_gcm, 1)
261      do i=1,nzy
262         co2y(i) = co2y(i) * nty(i)
263      enddo
264
265
266 
267
268c end                                           
269      return                                         
270      end 
271
272
273     
274     
275c***********************************************************************
276                                                           
277      subroutine NLTEdlvr09_CZALU(ig)
278
279c     jul 2011 malv+fgg
280c***********************************************************************
281                                                           
282      implicit none                                 
283                                                           
284!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! common variables and constants 
285                                                           
286      include 'nlte_paramdef.h'
287      include 'nlte_commons.h'
288     
289c arguments                           
290
291      integer  ig               !ADDED FOR TRACEBACK
292                                                           
293c local variables                               
294                                                           
295! matrixes and vectors                         
296                                                           
297      real*8 e110(nl), e210(nl), e310(nl), e410(nl) 
298      real*8 e121(nl), e112(nl)                     
299     
300      real*8 f1(nl,nl)
301                                                           
302      real*8 cax1(nl,nl), cax2(nl,nl), cax3(nl,nl)   
303      real*8 v1(nl), v2(nl), v3(nl)
304
305      real*8 alf11(nl,nl), alf12(nl,nl)             
306      real*8 alf21(nl,nl), alf31(nl,nl), alf41(nl,nl)           
307      real*8 a11(nl), a1112(nl,nl)                   
308      real*8            a1121(nl,nl), a1131(nl,nl), a1141(nl,nl)         
309      real*8 a21(nl), a2131(nl,nl), a2141(nl,nl)     
310      real*8            a2111(nl,nl), a2112(nl,nl)           
311      real*8 a31(nl), a3121(nl,nl), a3141(nl,nl)     
312      real*8            a3111(nl,nl), a3112(nl,nl)           
313      real*8 a41(nl), a4121(nl,nl), a4131(nl,nl)     
314      real*8            a4111(nl,nl), a4112(nl,nl)           
315      real*8 a12(nl), a1211(nl,nl)                   
316      real*8            a1221(nl,nl), a1231(nl,nl), a1241(nl,nl)         
317                                                           
318      real*8 aalf11(nl,nl),aalf21(nl,nl),aalf31(nl,nl),aalf41(nl,nl)         
319      real*8 aa11(nl), aa1121(nl,nl), aa1131(nl,nl), aa1141(nl,nl)           
320      real*8 aa21(nl), aa2111(nl,nl), aa2131(nl,nl), aa2141(nl,nl)           
321      real*8 aa31(nl), aa3111(nl,nl), aa3121(nl,nl), aa3141(nl,nl)           
322      real*8 aa41(nl), aa4111(nl,nl), aa4121(nl,nl), aa4131(nl,nl)           
323      real*8 aa12(nl)                             
324      real*8 aa1211(nl,nl), aa1221(nl,nl), aa1231(nl,nl), aa1241(nl,nl)     
325      real*8 aa1112(nl,nl), aa2112(nl,nl), aa3112(nl,nl), aa4112(nl,nl)     
326                                                           
327      real*8 aaalf11(nl,nl),aaalf21(nl,nl),aaalf31(nl,nl),
328     &     aaalf41(nl,nl)     
329      real*8 aaa11(nl),aaa1121(nl,nl),aaa1131(nl,nl),aaa1141(nl,nl)         
330      real*8 aaa21(nl),aaa2111(nl,nl),aaa2131(nl,nl),aaa2141(nl,nl)         
331      real*8 aaa31(nl),aaa3111(nl,nl),aaa3121(nl,nl),aaa3141(nl,nl)         
332      real*8 aaa41(nl),aaa4111(nl,nl),aaa4121(nl,nl),aaa4131(nl,nl)         
333                                                           
334      real*8 aaaalf11(nl,nl),aaaalf41(nl,nl)         
335      real*8 aaaa11(nl),aaaa1141(nl,nl)             
336      real*8 aaaa41(nl),aaaa4111(nl,nl)             
337                                                           
338                                                           
339                                                           
340! populations                                   
341      real*8 n10(nl), n11(nl)
342      real*8 n20(nl), n21(nl)
343      real*8 n30(nl), n31(nl)
344      real*8 n40(nl), n41(nl)
345     
346                                                           
347! productions and loses                         
348      real*8 d19a1,d19b1,d19c1
349      real*8 d19ap1,d19bp1,d19cp1 
350      real*8 d19a2,d19b2,d19c2
351      real*8 d19ap2,d19bp2,d19cp2 
352      real*8 d19a3,d19b3,d19c3
353      real*8 d19ap3,d19bp3,d19cp3 
354      real*8 d19a4,d19b4,d19c4
355      real*8 d19ap4,d19bp4,d19cp4 
356                                                           
357      real*8 l11, l12, l21, l31, l41                 
358      real*8 p11, p12, p21, p31, p41                 
359      real*8 p1112, p1211, p1221, p1231, p1241       
360      real*8 p1121, p1131, p1141                     
361      real*8 p2111, p2112, p2131, p2141             
362      real*8 p3111, p3112, p3121, p3141             
363      real*8 p4111, p4112, p4121, p4131             
364                                                           
365                                                           
366      real*8 ps11, ps21, ps31, ps41, ps12           
367     
368      real*8 pl11, pl12, pl21, pl31, pl41           
369                                                           
370c local constants and indexes                   
371                                                           
372      integer   ii              ! decides if output of tv,hr   
373      integer   icurt           ! decides if read/comp c.matrix
374
375      real*8 co2t                                   
376      real*8 ftest                                   
377                                                           
378      real*8 a11_einst(nl), a12_einst(nl)           
379      real*8 a21_einst(nl), a31_einst(nl), a41_einst(nl)         
380      real tsurf                                     
381
382      real*8 nu11, nu12, nu121, nu21, nu31, nu41
383                                                           
384      integer i, j, ik, isot , icurtishb                       
385      integer i_by15sh, i_col020, i_col010636                   
386                                                           
387                                                           
388c external functions and subroutines           
389                                                           
390      external planckdp                             
391      real*8    planckdp                               
392                                                           
393! subroutines called:                           
394!       mz4sub, dmzout, readc_mz4, mztf             
395                                                           
396!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!  start program
397                                                           
398
399      ii = 4
400      icurt = 1
401
402      call zero4v( aa11, aa21, aa31, aa41, nl)
403      call zero4m( aa1121, aa1131, aa1141, aalf11, nl)
404      call zero4m( aa2111, aa2131, aa2141, aalf21, nl)
405      call zero4m( aa3111, aa3121, aa3141, aalf31, nl)
406      call zero4m( aa4111, aa4121, aa4131, aalf41, nl)
407      call zero4m( aa1112, aa2112, aa3112, aa4112, nl)
408      call zero4m( aa1211, aa1221, aa1231, aa1241, nl)
409      call zero3v( aaa41, aaa31, aaa11, nl )
410      call zero3m( aaa4111, aaa4131, aaalf41, nl)
411      call zero3m( aaa3111, aaa3141, aaalf31, nl)
412      call zero3m( aaa1131, aaa1141, aaalf11, nl)
413      call zero2v( aaaa11, aaaa41, nl )
414      call zero2m( aaaa1141, aaaalf11, nl)
415      call zero2m( aaaa4111, aaaalf41, nl)
416     
417        !write (*,*)  ' --- c z a  simple ---    input_cza : ', input_cza   
418                                   
419
420      call zero3v (vt11,vt12,vt13,nl)               
421      call zero3v (vt21,vt22,vt23,nl)               
422      call zero3v (vt31,vt32,vt33,nl)               
423      call zero3v (vt41,vt42,vt43,nl)               
424                                                           
425      call zero3v (hr110,hr121,hr132,nl)             
426      call zero3v (hr210,hr221,hr232,nl)             
427      call zero3v (hr310,hr321,hr332,nl)             
428      call zero3v (hr410,hr421,hr432,nl)             
429      call zero3v (sl110,sl121,sl132,nl)             
430      call zero3v (sl210,sl221,sl232,nl)             
431      call zero3v (sl310,sl321,sl332,nl)             
432      call zero3v (sl410,sl421,sl432,nl)             
433     
434      call zero4v (el11,el21,el31,el41,nl)           
435      call zero4v (e110,e210,e310,e410,nl)           
436      call zero3v (el12,e121,e112,nl)               
437     
438      call zero3m (cax1,cax2,cax3,nl)               
439      call zerom (f1,nl)                             
440      call zero3v (v1,v2,v3,nl)                     
441     
442      call zero4m (alf11,alf21,alf31,alf41,nl)       
443      call zerom (alf12,nl)                         
444      call zero2v (a11,a12,nl)                       
445      call zero3v (a21,a31,a41,nl)                   
446     
447      call zero3m (a1121,a1131,a1141,nl)             
448      call zerom (a1112,nl)                         
449     
450      call zero3m (a1221,a1231,a1241,nl)             
451      call zerom (a1211,nl)                         
452     
453      call zero2m (a2111,a2112,nl)                   
454      call zero2m (a2131,a2141,nl)                   
455      call zero2m (a3111,a3112,nl)                   
456      call zero2m (a3121,a3141,nl)                   
457      call zero2m (a4111,a4112,nl)                   
458      call zero2m (a4121,a4131,nl)                   
459     
460                                                           
461      call zero4v (n11,n21,n31,n41,nl)               
462                                                           
463      nu11 = nu(1,1)                                 
464      nu12 = nu(1,2)                                 
465      nu121 = nu12-nu11                             
466     
467      nu21 = nu(2,1)                                 
468                                                           
469      nu31 = nu(3,1)                                 
470                                                           
471      nu41 = nu(4,1)                                 
472                                                           
473      ftest = 1.d0                                   
474      i_by15sh = 1
475      i_col020 = 1                                   
476                                                           
477      i_col010636 = 1                               
478                                                           
479                                                           
480 101  format(a1)                                 
481 180  format(a80)                                 
482                                                           
483                                                           
484c establishing molecular populations needed as input       
485      do i=1,nl                                     
486         n10(i) = dble( co2(i) * imr(1) )             
487         n20(i) = dble( co2(i) * imr(2) )             
488         n30(i) = dble( co2(i) * imr(3) )             
489         n40(i) = dble( co2(i) * imr(4) )             
490         if ( input_cza.ge.1 ) then                   
491            n11(i) = n10(i) *2.d0 *exp( dble(-ee*nu(1,1))/v626t1(i) )         
492            n21(i) = n20(i) *2.d0 *exp( dble(-ee*nu(2,1))/v628t1(i) )         
493            n31(i) = n30(i) *2.d0* exp( dble(-ee*nu(3,1))/v636t1(i) )         
494            n41(i) = n40(i) *2.d0* exp( dble(-ee*nu(4,1))/v627t1(i) )         
495         end if                                       
496      enddo                                   
497                                                       
498cc                                             
499cc   curtis matrix calculation                 
500cc                           
501      if ( input_cza.ge.1 ) then
502
503         if (itt_cza.eq.15 ) then
504           
505            call cm15um_hb_simple ( ig,icurt )
506           
507         elseif (itt_cza.eq.13) then
508           
509            call mztvc_626fh(ig)
510
511         endif
512
513      endif
514
515
516
517      do 4,i=nl,1,-1            !----------------------------------------------
518
519         co2t = dble ( co2(i) *(imr(1)+imr(3)+imr(2)+imr(4)) )     
520         
521         call getk ( t(i) )                             
522                                                                     
523         ps11 = 0.d0                                   
524         ps21 = 0.d0                                   
525         ps31 = 0.d0                                   
526         ps41 = 0.d0                                   
527         ps12 = 0.d0 
528                                 
529         ! V-T productions and losses V-T
530                                                           
531         isot = 1
532         d19b1 = dble(k19ba(isot)*co2t+k19bb(isot)*n2(i))         
533     @        + dble(k19bc(isot)*co(i))                   
534         d19c1 = dble(k19ca(isot)*co2t+k19cb(isot)*n2(i))         
535     @        + dble(k19cc(isot)*co(i))                   
536         d19bp1 = dble( k19bap(isot)*co2t + k19bbp(isot)*n2(i) ) 
537     @        + dble( k19bcp(isot)*co(i) )                 
538         d19cp1 = dble( k19cap(isot)*co2t + k19cbp(isot)*n2(i) ) 
539     @        + dble( k19ccp(isot)*co(i) )                 
540         isot = 2
541         d19c2 =dble(k19ca(isot)*co2t+k19cb(isot)*n2(i))         
542     @        + dble(k19cc(isot)*co(i))                   
543         d19cp2 =dble( k19cap(isot)*co2t + k19cbp(isot)*n2(i) )   
544     @        + dble( k19ccp(isot)*co(i) )                 
545         isot = 3
546         d19c3 =dble(k19ca(isot)*co2t+k19cb(isot)*n2(i))         
547     @        + dble(k19cc(isot)*co(i))                   
548         d19cp3 =dble( k19cap(isot)*co2t + k19cbp(isot)*n2(i) )   
549     @        + dble( k19ccp(isot)*co(i) )               
550         isot = 4
551         d19c4 =dble(k19ca(isot)*co2t+k19cb(isot)*n2(i))         
552     @        + dble(k19cc(isot)*co(i))                   
553         d19cp4 =dble( k19cap(isot)*co2t + k19cbp(isot)*n2(i) )   
554     @        + dble(k19ccp(isot)*co(i) )                 
555                                !
556         l11 = d19c1 + k20c(1)*dble(o3p(i))             
557         p11 = ( d19cp1 + k20cp(1)*dble(o3p(i)) ) * n10(i)         
558         l21 = d19c2 + k20c(2)*dble(o3p(i))             
559         p21 = ( d19cp2 + k20cp(2)*dble(o3p(i)) ) *n20(i)           
560         l31 = d19c3 + k20c(3)*dble(o3p(i))             
561         p31 = ( d19cp3 + k20cp(3)*dble(o3p(i)) ) *n30(i)           
562         l41 = d19c4 + k20c(4)*dble(o3p(i))             
563         p41 = ( d19cp4 + k20cp(4)*dble(o3p(i)) ) *n40(i)           
564           
565          ! Addition of V-V
566       
567         l11 = l11 + k21cp(2)*n20(i) + k21cp(3)*n30(i) + k21cp(4)*n40(i)     
568         p1121 = k21c(2) * n10(i)                     
569         p1131 = k21c(3) * n10(i)                     
570         p1141 = k21c(4) * n10(i)                     
571          !
572         l21 = l21 + k21c(2)*n10(i) + k23k21c*n30(i) + k24k21c*n40(i)         
573         p2111 = k21cp(2) * n20(i)                   
574         p2131 = k23k21cp * n20(i)                   
575         p2141 = k24k21cp * n20(i)                   
576          !                                                 
577         l31 = l31 + k21c(3)*n10(i) + k23k21cp*n20(i) + k34k21c*n40(i)       
578         p3111 = k21cp(3)* n30(i)                     
579         p3121 = k23k21c * n30(i)                     
580         p3141 = k34k21cp* n30(i)                     
581          !                                                 
582         l41 = l41 + k21c(4)*n10(i) + k24k21cp*n20(i) + k34k21cp*n30(i)       
583         p4111 = k21cp(4)* n40(i)                     
584         p4121 = k24k21c * n40(i)                     
585         p4131 = k34k21c * n40(i)                     
586                                                           
587                                                           
588         if ( input_cza.ge.1 ) then         
589                                                           
590            l12 = d19b1                                   
591     @           + k20b(1)*dble(o3p(i))                       
592     @           + k21b(1)*n10(i)                             
593     @           + k33c*( n20(i) + n30(i) + n40(i) )         
594            p12 = k21bp(1)*n11(i) * n11(i)               
595            p1211 = d19bp1 + k20bp(1)*dble(o3p(i))       
596            p1221 = k33cp(2)*n11(i)                       
597            p1231 = k33cp(3)*n11(i)                       
598            p1241 = k33cp(4)*n11(i)                       
599                                                           
600            l11 = l11 + d19bp1                           
601     @           + k20bp(1)*dble(o3p(i))                 
602     @           + 2.d0 * k21bp(1) * n11(i)               
603     @           +   k33cp(2)*n21(i) + k33cp(3)*n31(i) + k33cp(4)*n41(i)
604            p1112 = d19b1                               
605     @           + k20b(1)*dble(o3p(i))                   
606     @           + 2.d0*k21b(1)*n10(i)                   
607     @           + k33c*( n20(i) + n30(i) + n40(i) )     
608                                                           
609            l21 = l21 + k33cp(2)*n11(i)                   
610            p2112 = k33c*n20(i)                           
611           
612            l31 = l31 + k33cp(3)*n11(i)                   
613            p3112 = k33c*n30(i)                           
614                                                           
615            l41 = l41 + k33cp(4)*n11(i)                   
616            p4112 = k33c*n40(i)                           
617                                                           
618         end if                                         
619                                                           
620
621          ! Changes in local losses for ITT=13,15 cases
622
623         a21_einst(i) = 1.3452d00 * 1.8 / 4.0 * taustar21(i)   
624         a31_einst(i) = 1.1878d00 * 1.8 / 4.0 * taustar31(i)   
625         a41_einst(i) = 1.2455d00 * 1.8 / 4.0 * taustar41(i)   
626
627         l21 = l21 + a21_einst(i)             
628         l31 = l31 + a31_einst(i)             
629         l41 = l41 + a41_einst(i)             
630         
631         if (input_cza.ge.1 .and. itt_cza.eq.13) then
632            a12_einst(i) = 4.35d00 / 3.0d0 * 1.8 / 4.0 * taustar12(i)
633            l12=l12+a12_einst(i) 
634         endif
635
636         if (itt_cza.eq.24) then
637            a11_einst(i) = a11_einst(i)  * 1.8 / 4.0 * taustar11(i)
638            l11 = l11 + a11_einst(i)
639         endif
640           
641
642          !  vectors and matrices for the formulation                 
643
644         a11(i) = dble(gamma*nu11**3.) * 1.d0/2.d0 * (p11+ps11) /
645     @               (n10(i)*l11) 
646         a1121(i,i) = dble((nu11/nu21))**3.d0 * n20(i)/n10(i) *p1121/l11
647         a1131(i,i) = dble((nu11/nu31))**3.d0 * n30(i)/n10(i) *p1131/l11
648         a1141(i,i) = dble((nu11/nu41))**3.d0 * n40(i)/n10(i) *p1141/l11
649         e110(i) = 2.d0* dble(vlight*nu11**2.) * 1.d0/2.d0 /
650     @        ( n10(i) * l11 )   
651                                                           
652         a21(i) = dble( gamma*nu21**3.) * 1.d0/2.d0 *
653     @        (p21+ps21)/(n20(i)*l21)   
654         a2111(i,i) = dble((nu21/nu11))**3.d0 * n10(i)/n20(i) *p2111/l21
655         a2131(i,i) = dble((nu21/nu31))**3.d0 * n30(i)/n20(i) *p2131/l21   
656         a2141(i,i) = dble((nu21/nu41))**3.d0 * n40(i)/n20(i) *p2141/l21   
657         e210(i) = 2.d0*dble(vlight*nu21**2.) * 1.d0/2.d0 /
658     @        ( n20(i) * l21 )   
659                                                           
660         a31(i) = dble(gamma*nu31**3.) * 1.d0/2.d0 * (p31+ps31) /
661     @        (n30(i)*l31) 
662         a3111(i,i) = dble((nu31/nu11))**3.d0 * n10(i)/n30(i) *p3111/l31   
663         a3121(i,i) = dble((nu31/nu21))**3.d0 * n20(i)/n30(i) *p3121/l31   
664         a3141(i,i) = dble((nu31/nu41))**3.d0 * n40(i)/n30(i) *p3141/l31   
665         e310(i) = 2.d0*dble(vlight*nu31**2.) * 1.d0/2.d0 /
666     @        ( n30(i) * l31 )   
667         
668         a41(i) = dble(gamma*nu41**3.) * 1.d0/2.d0 * (p41+ps41) /
669     @        (n40(i)*l41) 
670         a4111(i,i) = dble((nu41/nu11))**3.d0 * n10(i)/n40(i) *p4111/l41   
671         a4121(i,i) = dble((nu41/nu21))**3.d0 * n20(i)/n40(i) *p4121/l41   
672         a4131(i,i) = dble((nu41/nu31))**3.d0 * n30(i)/n40(i) *p4131/l41
673         e410(i) = 2.d0*dble(vlight*nu41**2.) * 1.d0/2.d0 /
674     @        ( n40(i) * l41 )   
675                                                           
676         if (input_cza.ge.1) then                       
677           
678            a1112(i,i) = dble((nu11/nu121))**3.d0 * n11(i)/n10(i) *
679     @           p1112/l11   
680            a2112(i,i) = dble((nu21/nu121))**3.d0 * n11(i)/n20(i) *
681     @           p2112/l21   
682            a3112(i,i) = dble((nu31/nu121))**3.d0 * n11(i)/n30(i) *
683     @           p3112/l31   
684            a4112(i,i) = dble((nu41/nu121))**3.d0 * n11(i)/n40(i) *
685     @           p4112/l41   
686            e112(i) = -2.d0*dble(vlight*nu11**3.)/nu121 /2.d0 /
687     @           ( n10(i)*l11 )   
688            a12(i) = dble( gamma*nu121**3.) *2.d0/4.d0* (p12+ps12)/
689     @           (n11(i)*l12) 
690            a1211(i,i) = dble((nu121/nu11))**3.d0 * n10(i)/n11(i) *
691     @           p1211/l12   
692            a1221(i,i) = dble((nu121/nu21))**3.d0 * n20(i)/n11(i) *
693     @           p1221/l12   
694            a1231(i,i) = dble((nu121/nu31))**3.d0 * n30(i)/n11(i) *
695     @           p1231/l12   
696            a1241(i,i) = dble((nu121/nu41))**3.d0 * n40(i)/n11(i) *
697     @           p1241/l12   
698            e121(i) = 2.d0*dble(vlight*nu121**2.) *2.d0/4.d0 /
699     @           ( n11(i) * l12 ) 
700                                                           
701         end if                                         
702                                                           
703                                                           
704 4    continue    !-------------------------------------------------------   
705                                                           
706                                                           
707        ! Change C.M.
708                                                           
709      do i=1,nl                                   
710         do j=1,nl                                 
711            c210(i,j) = 0.0d0                             
712            c310(i,j) = 0.0d0                             
713            c410(i,j) = 0.0d0                             
714         end do                                     
715      end do
716      if ( itt_cza.eq.13 ) then
717         do i=1,nl                                   
718            do j=1,nl                                 
719               c121(i,j) = 0.0d0         
720            end do                                     
721         end do
722      endif
723        !Añadido para hacer diagonal C121
724!        if ( itt_cza.eq.15 ) then
725!         do i=1,nl                                   
726!           do j=1,nl
727!               if(abs(i-j).eq.1.or.abs(i-j).eq.2) c121(i,j) = 0.0d0         
728!           end do                                     
729!         end do
730!        endif
731      if ( itt_cza.eq.24 ) then
732         do i=1,nl                                   
733            do j=1,nl                                 
734               c110(i,j) = 0.0d0         
735            end do                                     
736         end do
737      endif
738                                                           
739        ! Lower Boundary
740      tsurf = t(1) + tsurf_excess                   
741      do i=1,nl                                     
742         sl110(i) = sl110(i) + vc110(i) * planckdp( tsurf, nu11 )
743         sl210(i) = sl210(i) + vc210(i) * planckdp( tsurf, nu21 )
744         sl310(i) = sl310(i) + vc310(i) * planckdp( tsurf, nu31 )
745         sl410(i) = sl410(i) + vc410(i) * planckdp( tsurf, nu41 )
746      end do                                         
747      if (input_cza.ge.1) then                       
748         do i=1,nl                                     
749            sl121(i) = sl121(i) + vc121(i) * planckdp( tsurf, nu121 )
750         end do     
751      endif
752               
753                                             
754        !!!!!!!!!!!! Solucion del sistema
755                                                           
756        !! Paso 0 :  Calculo de los alphas   alf11, alf21, alf31, alf41, alf12
757
758      call unit  ( cax2, nl )
759                                         
760      call diago ( cax1, e110, nl )
761      call mulmm ( cax3, cax1,c110, nl )                           
762!        cax3=matmul(cax1,c110)
763      call resmm ( alf11, cax2,cax3, nl )                           
764
765      call diago ( cax1, e210, nl )     
766      call mulmm ( cax3, cax1,c210, nl )                             
767!        cax3=matmul(cax1,c210)
768      call resmm ( alf21, cax2,cax3, nl )                           
769
770      call diago ( cax1, e310, nl )     
771      call mulmm ( cax3, cax1,c310, nl )                             
772!        cax3=matmul(cax1,c310)
773      call resmm ( alf31, cax2,cax3, nl )                           
774        !
775      call diago ( cax1, e410, nl )     
776      call mulmm ( cax3, cax1,c410, nl )                           
777!        cax3=matmul(cax1,c410)
778      call resmm ( alf41, cax2,cax3, nl )                         
779           !
780!        if(ig.eq.2223.and.input_cza.eq.1) then
781!           open(168,file='output_curtis_c121diagminus2.dat')
782!           do i=1,nl
783!              do j=1,nl
784!                 write(168,*)i,j,c110(i,j),c121(i,j)
785!              enddo
786!           enddo
787!           close(168)
788!           open(178,file='output_taustar.dat')
789!           do i=1,nl
790!              write(178,*)i,taustar21(i),taustar31(i),taustar41(i)
791!           enddo
792!           close(178)
793!        endif
794      if (input_cza.ge.1) then                   
795         call diago ( cax1, e121, nl ) 
796         call mulmm ( cax3, cax1,c121, nl )                       
797!        cax3=matmul(cax1,c121)
798         call resmm ( alf12, cax2,cax3, nl )
799      endif
800 
801        !! Paso 1 :  Calculo de vectores y matrices con 1 barra (aa***)
802           
803      if (input_cza.eq.0) then  !  Skip paso 1, pues el12 no se calcula
804
805          ! el11
806         call sypvvv( aa11, a11,e110,sl110, nl )
807         call samem( aa1121, a1121, nl )
808         call samem( aa1131, a1131, nl )
809         call samem( aa1141, a1141, nl )
810         call samem( aalf11, alf11, nl )
811         
812          ! el21
813         call sypvvv( aa21, a21,e210,sl210, nl )
814         call samem( aa2111, a2111, nl )
815         call samem( aa2131, a2131, nl )
816         call samem( aa2141, a2141, nl )
817         call samem( aalf21, alf21, nl )
818
819          ! el31
820         call sypvvv( aa31, a31,e310,sl310, nl )
821         call samem( aa3111, a3111, nl )
822         call samem( aa3121, a3121, nl )
823         call samem( aa3141, a3141, nl )
824         call samem( aalf31, alf31, nl )
825
826          ! el41
827         call sypvvv( aa41, a41,e410,sl410, nl )
828         call samem( aa4111, a4111, nl )
829         call samem( aa4121, a4121, nl )
830         call samem( aa4131, a4131, nl )
831         call samem( aalf41, alf41, nl )
832
833
834      else                      !      (input_cza.ge.1) ,   FH !
835
836
837         call sypvvv( v1, a12,e121,sl121, nl ) ! a12 + e121 * sl121
838
839          ! aa11
840         call sypvvv( v2, a11,e110,sl110, nl )
841         call trucommvv( aa11 , alf12,a1112,v2, v1, nl )
842           
843          ! aalf11
844         call invdiag( cax1, a1112, nl )
845         call mulmm( cax2, alf12, cax1, nl ) ! alf12 * (1/a1112)
846!          cax2=matmul(alf12,cax1)
847         call mulmm( cax3, cax2, alf11, nl )
848!          cax3=matmul(cax2,alf11)
849         
850         call resmm( aalf11, cax3, a1211, nl )
851          ! aa1121
852         call trucodiag(aa1121, alf12,a1112,a1121, a1221, nl)
853          ! aa1131
854         call trucodiag(aa1131, alf12,a1112,a1131, a1231, nl)
855          ! aa1141
856         call trucodiag(aa1141, alf12,a1112,a1141, a1241, nl)
857
858           
859          ! aa21
860         call sypvvv( v2, a21,e210,sl210, nl )
861         call trucommvv( aa21 , alf12,a2112,v2, v1, nl )
862
863          ! aalf21
864         call invdiag( cax1, a2112, nl )
865         call mulmm( cax2, alf12, cax1, nl ) ! alf12 * (1/a2112)
866!          cax2=matmul(alf12,cax1)
867         call mulmm( cax3, cax2, alf21, nl )
868!          cax3=matmul(cax2,alf21)
869         call resmm( aalf21, cax3, a1221, nl )
870          ! aa2111
871         call trucodiag(aa2111, alf12,a2112,a2111, a1211, nl)
872          ! aa2131
873         call trucodiag(aa2131, alf12,a2112,a2131, a1231, nl)
874          ! aa2141
875         call trucodiag(aa2141, alf12,a2112,a2141, a1241, nl)
876
877         
878          ! aa31
879         call sypvvv( v2, a31,e310,sl310, nl )
880         call trucommvv( aa31 , alf12,a3112,v2, v1, nl )
881          ! aalf31
882         call invdiag( cax1, a3112, nl )
883         call mulmm( cax2, alf12, cax1, nl ) ! alf12 * (1/a3112)
884!          cax2=matmul(alf12,cax1)
885         call mulmm( cax3, cax2, alf31, nl )
886!          cax3=matmul(cax2,alf31)
887         call resmm( aalf31, cax3, a1231, nl )
888          ! aa3111
889         call trucodiag(aa3111, alf12,a3112,a3111, a1211, nl)
890          ! aa3121
891         call trucodiag(aa3121, alf12,a3112,a3121, a1221, nl)
892          ! aa3141
893         call trucodiag(aa3141, alf12,a3112,a3141, a1241, nl)
894 
895
896          ! aa41
897         call sypvvv( v2, a41,e410,sl410, nl )
898         call trucommvv( aa41 , alf12,a4112,v2, v1, nl )
899          ! aalf41
900         call invdiag( cax1, a4112, nl )
901         call mulmm( cax2, alf12, cax1, nl ) ! alf12 * (1/a4112)
902!          cax2=matmul(alf12,cax1)
903         call mulmm( cax3, cax2, alf41, nl )
904!          cax3=matmul(cax2,alf41)
905         call resmm( aalf41, cax3, a1241, nl )
906          ! aa4111
907         call trucodiag(aa4111, alf12,a4112,a4111, a1211, nl)
908          ! aa4121
909         call trucodiag(aa4121, alf12,a4112,a4121, a1221, nl)
910          ! aa4131
911         call trucodiag(aa4131, alf12,a4112,a4131, a1231, nl)
912
913      endif                     ! Final  caso input_cza.ge.1
914
915
916         !! Paso 2 :  Calculo de vectores y matrices con 2 barras (aaa***)
917
918         ! aaalf41
919      call invdiag( cax1, aa4121, nl )
920      call mulmm( cax2, aalf21, cax1, nl ) ! alf21 * (1/a4121)
921!         cax2=matmul(aalf21,cax1)
922      call mulmm( cax3, cax2, aalf41, nl )
923!         cax3=matmul(cax2,aalf41)
924      call resmm( aaalf41, cax3, aa2141, nl )
925         ! aaa41
926      call trucommvv(aaa41, aalf21,aa4121,aa41, aa21, nl)
927         ! aaa4111
928      call trucodiag(aaa4111, aalf21,aa4121,aa4111, aa2111, nl)
929         ! aaa4131
930      call trucodiag(aaa4131, aalf21,aa4121,aa4131, aa2131, nl)
931
932         ! aaalf31
933      call invdiag( cax1, aa3121, nl )
934      call mulmm( cax2, aalf21, cax1, nl ) ! alf21 * (1/a3121)
935!         cax2=matmul(aalf21,cax1)
936      call mulmm( cax3, cax2, aalf31, nl )
937!         cax3=matmul(cax2,aalf31)
938      call resmm( aaalf31, cax3, aa2131, nl )
939         ! aaa31
940      call trucommvv(aaa31, aalf21,aa3121,aa31, aa21, nl)
941         ! aaa3111
942      call trucodiag(aaa3111, aalf21,aa3121,aa3111, aa2111, nl)
943         ! aaa3141
944      call trucodiag(aaa3141, aalf21,aa3121,aa3141, aa2141, nl)
945
946         ! aaalf11
947      call invdiag( cax1, aa1121, nl )
948      call mulmm( cax2, aalf21, cax1, nl ) ! alf21 * (1/a1121)
949!         cax2=matmul(aalf21,cax1)
950      call mulmm( cax3, cax2, aalf11, nl )
951!         cax3=matmul(cax2,aalf11)
952      call resmm( aaalf11, cax3, aa2111, nl )
953         ! aaa11
954      call trucommvv(aaa11, aalf21,aa1121,aa11, aa21, nl)
955         ! aaa1131
956      call trucodiag(aaa1131, aalf21,aa1121,aa1131, aa2131, nl)
957         ! aaa1141
958      call trucodiag(aaa1141, aalf21,aa1121,aa1141, aa2141, nl)
959
960
961         !! Paso 3 :  Calculo de vectores y matrices con 3 barras (aaaa***)
962
963         ! aaaalf41
964      call invdiag( cax1, aaa4131, nl )
965      call mulmm( cax2, aaalf31, cax1, nl ) ! aaalf31 * (1/aaa4131)
966!         cax2=matmul(aaalf31,cax1)
967      call mulmm( cax3, cax2, aaalf41, nl )
968!         cax3=matmul(cax2,aaalf41)
969      call resmm( aaaalf41, cax3, aaa3141, nl )
970         
971         ! aaaa41
972      call trucommvv(aaaa41, aaalf31,aaa4131,aaa41, aaa31, nl)
973         ! aaaa4111
974      call trucodiag(aaaa4111, aaalf31,aaa4131,aaa4111,aaa3111, nl)
975
976         ! aaaalf11
977      call invdiag( cax1, aaa1131, nl )
978      call mulmm( cax2, aaalf31, cax1, nl ) ! aaalf31 * (1/aaa4131)
979!         cax2=matmul(aaalf31,cax1)
980      call mulmm( cax3, cax2, aaalf11, nl )
981!     cax3=matmul(cax2,aaalf11)
982      call resmm( aaaalf11, cax3, aaa3111, nl )
983         ! aaaa11
984      call trucommvv(aaaa11, aaalf31,aaa1131,aaa11, aaa31, nl)
985         ! aaaa1141
986      call trucodiag(aaaa1141, aaalf31,aaa1131,aaa1141,aaa3141, nl)
987
988
989         !! Paso 4 :  Calculo de vectores y matrices finales y calculo de J1
990
991      call trucommvv(v1, aaaalf41,aaaa1141,aaaa11, aaaa41, nl)
992         !
993      call invdiag( cax1, aaaa1141, nl )
994      call mulmm( cax2, aaaalf41, cax1, nl ) ! aaaalf41 * (1/aaaa1141)
995!         cax2=matmul(aaaalf41,cax1)
996      call mulmm( cax3, cax2, aaaalf11, nl )
997!         cax3=matmul(cax2,aaaalf11)
998      call resmm( cax1, cax3, aaaa4111, nl )
999         !
1000      call LUdec ( el11, cax1, v1, nl, nl2 )
1001
1002         ! Solucion para el41
1003      call sypvmv( v1, aaaa41, aaaa4111,el11, nl )
1004      call LUdec ( el41, aaaalf41, v1, nl, nl2 )
1005
1006         ! Solucion para el31
1007      call sypvmv( v2, aaa31, aaa3111,el11, nl )
1008      call sypvmv( v1,    v2, aaa3141,el41, nl )
1009      call LUdec ( el31, aaalf31, v1, nl, nl2 )
1010
1011         ! Solucion para el21
1012      call sypvmv( v3, aa21, aa2111,el11, nl )
1013      call sypvmv( v2,   v3, aa2131,el31, nl )
1014      call sypvmv( v1,   v2, aa2141,el41, nl )
1015      call LUdec ( el21, aalf21, v1, nl, nl2 )
1016
1017         !!!
1018      el11(1) = planckdp( t(1), nu11 )         
1019      el21(1) = planckdp( t(1), nu21 )         
1020      el31(1) = planckdp( t(1), nu31 )         
1021      el41(1) = planckdp( t(1), nu41 )         
1022      el11(nl) = 2.d0 * el11(nl-1) - el11(nl2)   
1023      el21(nl) = 2.d0 * el21(nl-1) - el21(nl2)   
1024      el31(nl) = 2.d0 * el31(nl-1) - el31(nl2)   
1025      el41(nl) = 2.d0 * el41(nl-1) - el41(nl2)   
1026                                                           
1027      call mulmv ( v1, c110,el11, nl )               
1028      call sumvv ( hr110, v1,sl110, nl )             
1029
1030         ! Solucion para el12
1031      if (input_cza.ge.1) then   
1032
1033         call sypvmv( v1, a12, a1211,el11, nl )
1034         call sypvmv( v3,  v1, a1221,el21, nl )
1035         call sypvmv( v2,  v3, a1231,el31, nl )
1036         call sypvmv( v1,  v2, a1241,el41, nl )
1037         call LUdec ( el12, alf12, v1, nl, nl2 )
1038
1039         el12(1) = planckdp( t(1), nu121 )           
1040         el12(nl) = 2.d0 * el12(nl-1) - el12(nl2)   
1041
1042         if (itt_cza.eq.15) then
1043            call mulmv ( v1, c121,el12, nl )           
1044            call sumvv ( hr121, v1,sl121, nl )           
1045         endif
1046         
1047      end if                                       
1048                                                           
1049                                                           
1050                                                           
1051      if (input_cza.lt.1) then
1052
1053         do i=1,nl                                                           
1054            pl11 = el11(i)/dble( gamma * nu11**3.0d0  * 1./2. / n10(i) )   
1055            pl21 = el21(i)/dble( gamma * nu21**3.0d0  * 1./2. / n20(i) )   
1056            pl31 = el31(i)/dble( gamma * nu31**3.0d0  * 1./2. / n30(i) )   
1057            pl41 = el41(i)/dble( gamma * nu41**3.0d0  * 1./2. / n40(i) )   
1058            vt11(i) = dble(-ee*nu11) / log( abs(pl11) / (2.0d0*n10(i)) )   
1059            vt21(i) = dble(-ee*nu21) / log( abs(pl21) / (2.0d0*n20(i)) )   
1060            vt31(i) = dble(-ee*nu31) / log( abs(pl31) / (2.0d0*n30(i)) )   
1061            vt41(i) = dble(-ee*nu41) / log( abs(pl41) / (2.0d0*n40(i)) )
1062            hr210(i) = sl210(i) - hplanck*vlight*nu21 * a21_einst(i)*pl21
1063            hr310(i) = sl310(i) - hplanck*vlight*nu31 * a31_einst(i)*pl31
1064            hr410(i) = sl410(i) - hplanck*vlight*nu41 * a41_einst(i)*pl41
1065!            hr410(i) = 0.
1066         enddo
1067
1068         call dinterconnection ( v626t1, vt11 )         
1069         call dinterconnection ( v628t1, vt21 )         
1070         call dinterconnection ( v636t1, vt31 )         
1071         call dinterconnection ( v627t1, vt41 )         
1072
1073      else
1074                                               
1075         do i=1,nl                                                           
1076            pl21 = el21(i)/dble( gamma * nu21**3.0d0  * 1./2. / n20(i) )   
1077            pl31 = el31(i)/dble( gamma * nu31**3.0d0  * 1./2. / n30(i) )   
1078            pl41 = el41(i)/dble( gamma * nu41**3.0d0  * 1./2. / n40(i) )
1079            hr210(i) = sl210(i) - hplanck*vlight*nu21 * a21_einst(i)*pl21
1080            hr310(i) = sl310(i) - hplanck*vlight*nu31 * a31_einst(i)*pl31
1081            hr410(i) = sl410(i) - hplanck*vlight*nu41 * a41_einst(i)*pl41
1082!            hr410(i) = 0.
1083            if (itt_cza.eq.13) then                   
1084               pl12 = el12(i)/dble(gamma*nu121**3.0d0*2./4./n11(i)) 
1085               hr121(i) = - hplanck*vlight * nu121 * a12_einst(i) * pl12       
1086               hr121(i) = hr121(i) + sl121(i)
1087            endif                         
1088         enddo
1089
1090      endif
1091
1092        ! K/Dday
1093      do i=1,nl                                     
1094         hr110(i)=hr110(i)*( hrkday_factor(i) / nt(i) )
1095         hr210(i)=hr210(i)*( hrkday_factor(i) / nt(i) )           
1096         hr310(i)=hr310(i)*( hrkday_factor(i) / nt(i) )           
1097         hr410(i)=hr410(i)*( hrkday_factor(i) / nt(i) )           
1098         hr121(i)=hr121(i)*( hrkday_factor(i) / nt(i) )           
1099      end do                                         
1100                                                           
1101                                                           
1102
1103c  output                                       
1104                                                           
1105        !codigo = codeout                                                     
1106        !call dmzout_tv ( 1 )             
1107        !call dmzout_hr ( 1 )             
1108
1109c final subrutina                                                           
1110      return                                         
1111      end   
1112
1113c***********************************************************************
1114c       hrkday_convert.f                             
1115c                                               
1116c       fortran function that returns the factor for conversion from         
1117c       hr' [erg s-1 cm-3] to hr [ k day-1 ]           
1118c
1119c       mar 2010        fgg      adapted to GCM
1120c       jan 99          malv     add o2 as major component.
1121c       ago 98          malv     also returns cp_avg,pm_avg
1122c       jul 98          malv     first version.                 
1123c***********************************************************************
1124                                               
1125        function hrkday_convert                       
1126     @     ( mmean_nlte,cpmean_nlte )         
1127                                               
1128        implicit none                           
1129                         
1130        include 'comcstfi.h'
1131        include 'param.h'
1132                                               
1133c argumentos                                   
1134        real    mmean_nlte,cpmean_nlte
1135        real    hrkday_convert                           
1136                                               
1137ccccccccccccccccccccccccccccccccccccc           
1138       
1139        hrkday_convert = daysec * n_avog /
1140     &                  ( cpmean_nlte * 1.e4 * mmean_nlte )
1141                                               
1142c end                                           
1143        return                                 
1144        end                                     
1145
1146c***********************************************************************
1147        subroutine sypvvv(a,b,c,d,n)
1148c       a(i)=b(i)+c(i)*d(i)
1149c       jul 2011 malv+fgg
1150c***********************************************************************
1151        real*8 a(n),b(n),c(n),d(n)
1152        integer n,i
1153        do 1,i=2,n-1
1154          a(i)= b(i) + c(i) * d(i)
1155 1      continue
1156        a(1) = 0.0d0
1157        a(n) = 0.0d0
1158        return
1159        end
1160
1161c***********************************************************************
1162        subroutine sypvmv(v,u,c,w,n)
1163c       inputs: matriz diagonal c , vectores u,w
1164c       output: vector v
1165c       Operacion a realizar:  v = u + c * w
1166
1167c       jul 2011 malv+fgg
1168c***********************************************************************
1169        real*8 v(n),u(n),c(n,n),w(n)
1170        integer n,i
1171        do 1,i=2,n-1
1172          v(i)= u(i) + c(i,i) * w(i)
1173 1      continue
1174        v(1) = 0.0d0
1175        v(n) = 0.0d0
1176        return
1177        end
1178
1179c***********************************************************************
1180        subroutine trucommvv(v,b,c,u,w,n)
1181c       inputs: matrices b,c , vectores u,w
1182c       output: vector v
1183c       Operacion a realizar:  v = b * c^(-1) * u + w
1184c       La matriz c va a ser invertida
1185c       c es diagonal, b no
1186c       Aprovechamos esa condicion para invertir c, y acelerar el calculo
1187c       jul 2011 malv+fgg 
1188c***********************************************************************
1189        real*8 v(n),b(n,n),c(n,n),u(n),w(n), sum
1190        integer n,i,j,k
1191        do 1,i=2,n-1
1192          sum=0.0d0
1193          do 2,j=2,n-1
1194            sum=sum+ (b(i,j)) * (u(j)/c(j,j))
1195 2        continue
1196          v(i) = sum + w(i)
1197 1      continue
1198        v(1) = 0.d0
1199        v(n) = 0.d0
1200        return
1201        end
1202
1203c***********************************************************************
1204        subroutine trucodiag(a,b,c,d,e,n)
1205c       inputs: matrices b,c,d,e
1206c       output: matriz diagonal a
1207c       Operacion a realizar:  a = b * c^(-1) * d + e
1208c       La matriz c va a ser invertida
1209c       Todas las matrices de entrada son diagonales excepto b
1210c       Aprovechamos esa condicion para invertir c, acelerar el calculo, y
1211c       ademas, para forzar que a sea diagonal
1212c       jul 2011 malv+fgg
1213c***********************************************************************
1214        real*8 a(n,n),b(n,n),c(n,n),d(n,n),e(n,n), sum
1215        integer n,i,j,k
1216        do 1,i=2,n-1
1217          sum=0.0d0
1218          do 2,j=2,n-1
1219            sum=sum+ (b(i,j)) * (d(j,j)/c(j,j))
1220 2        continue
1221          a(i,i) = sum + e(i,i)
1222 1      continue
1223        do k=1,n
1224          a(n,k) = 0.0d0
1225          a(1,k) = 0.0d0
1226          a(k,1) = 0.0d0
1227          a(k,n) = 0.0d0
1228        end do
1229        return
1230        end
1231
1232c***********************************************************************
1233        subroutine invdiag(a,b,n)
1234c       inverse of a diagonal matrix
1235c       jul 2011 malv
1236c***********************************************************************
1237        implicit none
1238
1239        integer n,i,j,k
1240        real*8 a(n,n),b(n,n)
1241
1242        do 1,i=2,n-1
1243          do 2,j=2,n-1
1244            if (i.eq.j) then
1245              a(i,j) = 1.d0/b(i,i)
1246            else
1247              a(i,j)=0.0d0
1248            end if
1249 2        continue
1250 1      continue
1251        do k=1,n
1252          a(n,k) = 0.0d0
1253          a(1,k) = 0.0d0
1254          a(k,1) = 0.0d0
1255          a(k,n) = 0.0d0
1256        end do
1257        return
1258        end
Note: See TracBrowser for help on using the repository browser.