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

Last change on this file since 524 was 498, checked in by emillour, 13 years ago

Mars GCM: Cleanup of the NLTE routines which have been packed together to limit the number of files.
Also enforced that file names are non-capitalized (needed by the create_make_gcm scripts to better evaluate dependencies when building the makefile).
FGG+EM

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