This file is indexed.

/usr/share/code_saturne/user_examples/usatch_4spe5reac.f90 is in code-saturne-data 4.2.0+repack-1build1.

This file is owned by root:root, with mode 0o644.

The actual contents of the file can be viewed below.

  1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
!-------------------------------------------------------------------------------

! This file is part of Code_Saturne, a general-purpose CFD tool.
!
! Copyright (C) 1998-2015 EDF S.A.
!
! This program is free software; you can redistribute it and/or modify it under
! the terms of the GNU General Public License as published by the Free Software
! Foundation; either version 2 of the License, or (at your option) any later
! version.
!
! This program is distributed in the hope that it will be useful, but WITHOUT
! ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
! FOR A PARTICULAR PURPOSE.  See the GNU General Public License for more
! details.
!
! You should have received a copy of the GNU General Public License along with
! this program; if not, write to the Free Software Foundation, Inc., 51 Franklin
! Street, Fifth Floor, Boston, MA 02110-1301, USA.

!-------------------------------------------------------------------------------

!> \file usatch_4spe5reac.f90
!> \brief Routines for user defined atmospheric chemical scheme
!>  Example of the computation of an atmopsheric chemical scheme assuming
!>  quasi steady equilibrium NOx scheme with 4 species and 5 reactions
!>  (equivalent to ichemistry = 1)
!> \remarks
!>  These routines should be generated by SPACK
!>  See CEREA: http://cerea.enpc.fr/polyphemus


!> kinetic
!> \brief Computation of kinetic rates for atmospheric chemistry
!------------------------------------------------------------------------------

!------------------------------------------------------------------------------
! Arguments
!------------------------------------------------------------------------------
!   mode          name              role
!------------------------------------------------------------------------------
!> \param[in]     nr                total number of chemical reactions
!> \param[in]     option_photolysis flag to activate or not photolysis reactions
!> \param[in]     azi               solar zenith angle
!> \param[in]     att               atmospheric attenuation variable
!> \param[in]     temp              temperature
!> \param[in]     press             pressure
!> \param[in]     xlw               water massic fraction
!> \param[out]    rk(nr)            kinetic rates
!______________________________________________________________________________

subroutine kinetic(nr,rk,temp,xlw,press,azi,att,                  &
     option_photolysis)

use entsor

implicit none

! Arguments

integer nr
double precision rk(nr),temp,xlw,press
double precision azi, att
integer option_photolysis

! Local variables

double precision effko,rapk,facteur,summ
double precision ylh2o

if(1.eq.1) then
  write(nfecra,9000)
  call csexit (1)
endif

! Example

!     Compute third body.
!     Conversion = Avogadro*1d-6/Perfect gas constant.
!     PRESS in Pascal, SUMM in molecules/cm3, TEMP in Kelvin

summ = press * 7.243d16 / temp

!     Number of water molecules computed from the massic fraction
!     (absolute humidity)

ylh2o = 29.d0*summ*xlw/(18.d0+11.d0*xlw)

!     For the zenithal angle at tropics

azi=abs(azi)

rk(  1) =  dexp(-0.8860689615829534d+02                           &
 - ( -0.5300000000000000d+03 )/temp)
rk(  1) = rk(  1) * summ * 0.2d0
rk(  2) =  dexp(-0.2653240882726044d+02                           &
 - (  0.1500000000000000d+04 )/temp)

if(option_photolysis.eq.2) then
 rk(  3)= 0.d0
elseif(option_photolysis.eq.1) then
if(azi.lt.0.d0)then
 stop
elseif(azi.ge. 0.00d+00 .and. azi.lt. 0.10d+02) then
 rk(  3)=-0.1302720567168795d-07
 rk(  3)=-0.7822279432831311d-06+(azi- 0.00d+00) * rk(  3)
 rk(  3)= 0.0000000000000000d+00+(azi- 0.00d+00) * rk(  3)
 rk(  3)= 0.9310260000000001d-02+(azi- 0.00d+00) * rk(  3)
elseif(azi.ge. 0.10d+02 .and. azi.lt. 0.20d+02) then
 rk(  3)= 0.3771617015067078d-08
 rk(  3)=-0.1173044113433769d-05+(azi- 0.10d+02) * rk(  3)
 rk(  3)=-0.1955272056716901d-04+(azi- 0.10d+02) * rk(  3)
 rk(  3)= 0.9219010000000000d-02+(azi- 0.10d+02) * rk(  3)
elseif(azi.ge. 0.20d+02 .and. azi.lt. 0.30d+02) then
 rk(  3)=-0.5859262388581815d-08
 rk(  3)=-0.1059895602981758d-05+(azi- 0.20d+02) * rk(  3)
 rk(  3)=-0.4188211773132428d-04+(azi- 0.20d+02) * rk(  3)
 rk(  3)= 0.8909950000000000d-02+(azi- 0.20d+02) * rk(  3)
elseif(azi.ge. 0.30d+02 .and. azi.lt. 0.40d+02) then
 rk(  3)=-0.7024567460738029d-08
 rk(  3)=-0.1235673474639213d-05+(azi- 0.30d+02) * rk(  3)
 rk(  3)=-0.6483780850753392d-04+(azi- 0.30d+02) * rk(  3)
 rk(  3)= 0.8379279999999999d-02+(azi- 0.30d+02) * rk(  3)
elseif(azi.ge. 0.40d+02 .and. azi.lt. 0.50d+02) then
 rk(  3)=-0.9202467768466835d-08
 rk(  3)=-0.1446410498461367d-05+(azi- 0.40d+02) * rk(  3)
 rk(  3)=-0.9165864823853972d-04+(azi- 0.40d+02) * rk(  3)
 rk(  3)= 0.7600310000000000d-02+(azi- 0.40d+02) * rk(  3)
elseif(azi.ge. 0.50d+02 .and. azi.lt. 0.60d+02) then
 rk(  3)=-0.1612556146540100d-07
 rk(  3)=-0.1722484531515342d-05+(azi- 0.50d+02) * rk(  3)
 rk(  3)=-0.1233475985383066d-03+(azi- 0.50d+02) * rk(  3)
 rk(  3)= 0.6529880000000000d-02+(azi- 0.50d+02) * rk(  3)
elseif(azi.ge. 0.60d+02 .and. azi.lt. 0.70d+02) then
 rk(  3)= 0.3226471363007382d-07
 rk(  3)=-0.2206251375477548d-05+(azi- 0.60d+02) * rk(  3)
 rk(  3)=-0.1626349576082332d-03+(azi- 0.60d+02) * rk(  3)
 rk(  3)= 0.5108030000000000d-02+(azi- 0.60d+02) * rk(  3)
elseif(azi.ge. 0.70d+02 .and. azi.lt. 0.78d+02) then
 rk(  3)= 0.2027078243961372d-06
 rk(  3)=-0.1238309966574737d-05+(azi- 0.70d+02) * rk(  3)
 rk(  3)=-0.1970805710287543d-03+(azi- 0.70d+02) * rk(  3)
 rk(  3)= 0.3293320000000000d-02+(azi- 0.70d+02) * rk(  3)
elseif(azi.ge. 0.78d+02 .and. azi.lt. 0.86d+02) then
 rk(  3)=-0.7448311471194499d-07
 rk(  3)= 0.3626677818932555d-05+(azi- 0.78d+02) * rk(  3)
 rk(  3)=-0.1779736282099126d-03+(azi- 0.78d+02) * rk(  3)
 rk(  3)= 0.1741210000000000d-02+(azi- 0.78d+02) * rk(  3)
elseif(azi.ge. 0.86d+02 .and. azi.lt. 0.90d+02) then
 rk(  3)= 0.2490309929270573d-05
 rk(  3)= 0.1839083065842406d-05+(azi- 0.86d+02) * rk(  3)
 rk(  3)=-0.1342475411316713d-03+(azi- 0.86d+02) * rk(  3)
 rk(  3)= 0.5113930000000000d-03+(azi- 0.86d+02) * rk(  3)
elseif(azi.ge.90.d0)then
 rk(  3)= 0.1632080000000000d-03
endif
if(att.lt.0.99999) rk(  3) = rk(  3) * att
endif

rk(  4) = summ * 6.0d-34 * (temp/3.d2) ** (-2.4d0)
rk(  4) = rk(  4) * summ * 0.2d0
rk(  5) =  dexp(-0.2590825451818744d+02                           &
 - ( -0.1800000000000000d+03 )/temp)

return

!--------
! FORMATS
!--------
#if defined(_CS_LANG_FR)
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : CHIMIE ATMOSPHERIQUE DEMANDEE'                   ,/,&
  '@    ========'                                                ,/,&
  '@    L''utilisateur a choisi de fournir son propre schema'    ,/,&
  '@    chimique. Or aucun schema n''a ete trouve.'              ,/,&
  '@'                                                            ,/,&
  '@  Le calcul ne sera pas execute.'                            ,/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#else
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : ATMOSPHERIC CHEMISTRY'                           ,/,&
  '@    ========'                                                ,/,&
  '@    The user choose to use its own chemical scheme'          ,/,&
  '@    However no scheme has been found'                        ,/,&
  '@'                                                            ,/,&
  '@  The computation will not be run                           ',/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#endif

end subroutine kinetic


!===============================================================================

!> fexchem
!> \brief Computation of the chemical production terms
!------------------------------------------------------------------------------

!------------------------------------------------------------------------------
! Arguments
!------------------------------------------------------------------------------
!   mode          name              role
!------------------------------------------------------------------------------
!> \param[in]     nr                total number of chemical reactions
!> \param[in]     ns                total number of chemical species
!> \param[in]     y                 concentrations vector
!> \param[in]     rk                kinetic rates
!> \param[in]     zcsourc           source term
!> \param[in]     convers_factor    conversion factors
!> \param[out]    chem              chemical production terms for every species
!______________________________________________________________________________

subroutine fexchem(ns,nr,y,rk,zcsourc,convers_factor,chem)

use entsor

implicit none

! Arguments

integer nr,ns
double precision rk(nr),y(ns),chem(ns),zcsourc(ns)
double precision convers_factor(ns)

! Local variables

integer i
double precision w(nr)
double precision conc(ns)

if (1.eq.1) then
  write(nfecra,9000)
  call csexit (1)
endif

! Example

do i=1,ns
  chem(i)=0.d0
enddo

!     Conversion mug/m3 to molecules/cm3.

do i = 1, ns
  conc(i) = y(i) * convers_factor(i)
enddo

!     Compute reaction rates.

call rates(ns,nr,rk,conc,w)

!     Chemical production terms.

chem(  3) = chem(  3) +  0.2000000000000000d+01 * w(  1)
chem(  4) = chem(  4) -  0.2000000000000000d+01 * w(  1)
chem(  2) = chem(  2) - w(  2)
chem(  3) = chem(  3) + w(  2)
chem(  4) = chem(  4) - w(  2)
chem(  1) = chem(  1) + w(  3)
chem(  3) = chem(  3) - w(  3)
chem(  4) = chem(  4) + w(  3)
chem(  1) = chem(  1) - w(  4)
chem(  2) = chem(  2) + w(  4)
chem(  1) = chem(  1) - w(  5)
chem(  3) = chem(  3) - w(  5)
chem(  4) = chem(  4) + w(  5)

!    Conversion molecules/cm3 to mug/m3.

do i = 1, ns
  chem(i) = chem(i) / convers_factor(i)
enddo

!     Volumic source terms.

do i=1,ns
  chem(i)=chem(i)+zcsourc(i)
enddo

return

!--------
! FORMATS
!--------
#if defined(_CS_LANG_FR)
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : CHIMIE ATMOSPHERIQUE DEMANDEE'                   ,/,&
  '@    ========'                                                ,/,&
  '@    L''utilisateur a choisi de fournir son propre schema'    ,/,&
  '@    chimique. Or aucun schema n''a ete trouve.'              ,/,&
  '@'                                                            ,/,&
  '@  Le calcul ne sera pas execute.'                            ,/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#else
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : ATMOSPHERIC CHEMISTRY'                           ,/,&
  '@    ========'                                                ,/,&
  '@    The user choose to use its own chemical scheme'          ,/,&
  '@    However no scheme has been found'                        ,/,&
  '@'                                                            ,/,&
  '@  The computation will not be run                           ',/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#endif

end subroutine fexchem


!===============================================================================

!> jacdchemdc
!> \brief Computation of the Jacobian matrix for atmospheric chemistry
!------------------------------------------------------------------------------

!------------------------------------------------------------------------------
! Arguments
!------------------------------------------------------------------------------
!   mode          name               role
!------------------------------------------------------------------------------
!> \param[in]     nr                 total number of chemical reactions
!> \param[in]     ns                 total number of chemical species
!> \param[in]     y                  concentrations vector
!> \param[in]     convers_factor     conversion factors of mug/m3 to
!>                                   molecules/cm3
!> \param[in]     convers_factor_jac conversion factors for the Jacobian matrix
!>                                   (Wmol(i)/Wmol(j))
!> \param[in]     rk                 kinetic rates
!> \param[out]    jacc               Jacobian matrix
!______________________________________________________________________________

subroutine jacdchemdc(ns,nr,y,convers_factor,                     &
                    convers_factor_jac,rk,jacc)

use entsor

implicit none

! Arguments

integer nr,ns
double precision rk(nr),y(ns),jacc(ns,ns)
double precision convers_factor(ns)
double precision convers_factor_jac(ns,ns)

! Local variables

integer i,j
double precision dw(nr,ns)
double precision conc(ns)

if(1.eq.1) then
  write(nfecra,9000)
  call csexit (1)
endif

! Example

do j=1,ns
 do i=1,ns
  jacc(i,j)=0.d0
 enddo
enddo

!     Conversion mug/m3 to molecules/cm3.

do i = 1, ns
   conc(i) = y(i) * convers_factor(i)
enddo

call dratedc(ns,nr,rk,conc,dw)

jacc(  3,  4) = jacc(  3,  4)+ 0.2000000000000000d+01*dw(  1,  4)
jacc(  3,  4) = jacc(  3,  4)+ 0.2000000000000000d+01*dw(  1,  4)
jacc(  4,  4) = jacc(  4,  4)- 0.2000000000000000d+01*dw(  1,  4)
jacc(  4,  4) = jacc(  4,  4)- 0.2000000000000000d+01*dw(  1,  4)
jacc(  2,  2) = jacc(  2,  2) - dw(  2,  2)
jacc(  2,  4) = jacc(  2,  4) - dw(  2,  4)
jacc(  3,  2) = jacc(  3,  2) + dw(  2,  2)
jacc(  3,  4) = jacc(  3,  4) + dw(  2,  4)
jacc(  4,  2) = jacc(  4,  2) - dw(  2,  2)
jacc(  4,  4) = jacc(  4,  4) - dw(  2,  4)
jacc(  1,  3) = jacc(  1,  3) + dw(  3,  3)
jacc(  3,  3) = jacc(  3,  3) - dw(  3,  3)
jacc(  4,  3) = jacc(  4,  3) + dw(  3,  3)
jacc(  1,  1) = jacc(  1,  1) - dw(  4,  1)
jacc(  2,  1) = jacc(  2,  1) + dw(  4,  1)
jacc(  1,  1) = jacc(  1,  1) - dw(  5,  1)
jacc(  1,  3) = jacc(  1,  3) - dw(  5,  3)
jacc(  3,  1) = jacc(  3,  1) - dw(  5,  1)
jacc(  3,  3) = jacc(  3,  3) - dw(  5,  3)
jacc(  4,  1) = jacc(  4,  1) + dw(  5,  1)
jacc(  4,  3) = jacc(  4,  3) + dw(  5,  3)

do j = 1, ns
  do i = 1, ns
    jacc(i,j) = jacc(i,j) * convers_factor_jac(i,j)
  enddo
enddo


return

!--------
! FORMATS
!--------
#if defined(_CS_LANG_FR)
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : CHIMIE ATMOSPHERIQUE DEMANDEE'                   ,/,&
  '@    ========'                                                ,/,&
  '@    L''utilisateur a choisi de fournir son propre schema'    ,/,&
  '@    chimique. Or aucun schema n''a ete trouve.'              ,/,&
  '@'                                                            ,/,&
  '@  Le calcul ne sera pas execute.'                            ,/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#else
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : ATMOSPHERIC CHEMISTRY'                           ,/,&
  '@    ========'                                                ,/,&
  '@    The user choose to use its own chemical scheme'          ,/,&
  '@    However no scheme has been found'                        ,/,&
  '@'                                                            ,/,&
  '@  The computation will not be run                           ',/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#endif

end subroutine jacdchemdc


!===============================================================================

!> rates
!> \brief Computation of reaction rates
!------------------------------------------------------------------------------

!------------------------------------------------------------------------------
! Arguments
!------------------------------------------------------------------------------
!   mode          name               role
!------------------------------------------------------------------------------
!> \param[in]     nr                 total number of chemical reactions
!> \param[in]     ns                 total number of chemical species
!> \param[in]     rk                 kinetic rates
!> \param[in]     y                  concentrations vector
!> \param[out]    w                  reaction rates
!______________________________________________________________________________

subroutine rates(ns,nr,rk,y,w)

use entsor

implicit none

! Arguments

integer nr,ns
double precision rk(nr),y(ns)
double precision w(nr)

! Local variables

if(1.eq.1) then
  write(nfecra,9000)
  call csexit (1)
endif

! Example

w(  1) =  rk(  1) * y(  4) * y(  4)
w(  2) =  rk(  2) * y(  2) * y(  4)
w(  3) =  rk(  3) * y(  3)
w(  4) =  rk(  4) * y(  1)
w(  5) =  rk(  5) * y(  1) * y(  3)

return

!--------
! FORMATS
!--------
#if defined(_CS_LANG_FR)
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : CHIMIE ATMOSPHERIQUE DEMANDEE'                   ,/,&
  '@    ========'                                                ,/,&
  '@    L''utilisateur a choisi de fournir son propre schema'    ,/,&
  '@    chimique. Or aucun schema n''a ete trouve.'              ,/,&
  '@'                                                            ,/,&
  '@  Le calcul ne sera pas execute.'                            ,/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#else
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : ATMOSPHERIC CHEMISTRY'                           ,/,&
  '@    ========'                                                ,/,&
  '@    The user choose to use its own chemical scheme'          ,/,&
  '@    However no scheme has been found'                        ,/,&
  '@'                                                            ,/,&
  '@  The computation will not be run                           ',/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#endif

end subroutine rates


!===============================================================================

!> dratedc
!> \brief Computation of derivatives of reaction rates
!------------------------------------------------------------------------------

!------------------------------------------------------------------------------
! Arguments
!------------------------------------------------------------------------------
!   mode          name               role
!------------------------------------------------------------------------------
!> \param[in]     nr                 total number of chemical reactions
!> \param[in]     ns                 total number of chemical species
!> \param[in]     rk                 kinetic rates
!> \param[in]     y                  concentrations vector
!> \param[out]    dw                 derivatives of reaction rates
!______________________________________________________________________________

subroutine dratedc(ns,nr,rk,y,dw)

use entsor

implicit none

! Arguments

integer nr,ns
double precision rk(nr),y(ns)
double precision dw(nr,ns)

! Local variables


if(1.eq.1) then
  write(nfecra,9000)
  call csexit (1)
endif

! Example

dw(  1,  4) =  rk(  1) * y(  4)
dw(  1,  4) =  rk(  1) * y(  4)
dw(  2,  2) =  rk(  2) * y(  4)
dw(  2,  4) =  rk(  2) * y(  2)
dw(  3,  3) =  rk(  3)
dw(  4,  1) =  rk(  4)
dw(  5,  1) =  rk(  5) * y(  3)
dw(  5,  3) =  rk(  5) * y(  1)

return

!--------
! FORMATS
!--------
#if defined(_CS_LANG_FR)
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : CHIMIE ATMOSPHERIQUE DEMANDEE'                   ,/,&
  '@    ========'                                                ,/,&
  '@    L''utilisateur a choisi de fournir son propre schema'    ,/,&
  '@    chimique. Or aucun schema n''a ete trouve.'              ,/,&
  '@'                                                            ,/,&
  '@  Le calcul ne sera pas execute.'                            ,/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#else
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : ATMOSPHERIC CHEMISTRY'                           ,/,&
  '@    ========'                                                ,/,&
  '@    The user choose to use its own chemical scheme'          ,/,&
  '@    However no scheme has been found'                        ,/,&
  '@'                                                            ,/,&
  '@  The computation will not be run                           ',/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#endif


end subroutine dratedc


!===============================================================================

!> lu_decompose
!> \brief Computation of LU factorization of matrix m
!------------------------------------------------------------------------------

!------------------------------------------------------------------------------
! Arguments
!------------------------------------------------------------------------------
!   mode          name               role
!------------------------------------------------------------------------------
!> \param[in]     ns                 matrix row number from the chemical
!>                                   species number
!> \param[in,out] m                  on entry, an invertible matrix.
!>                                   On exit, an LU factorization of m
!______________________________________________________________________________

subroutine lu_decompose (ns,m)

use entsor

implicit none

! Arguments

integer ns
double precision m(ns,ns)

! Local variables

double precision temp

if(1.eq.1) then
  write(nfecra,9000)
  call csexit (1)
endif

! Example

!     Upper part.
m(1, 3) = m(1, 3) / m(1, 1)

!     Upper part.
temp = m(2, 1) * m(1, 3)
m(2, 3) = ( m(2, 3) - temp ) / m(2, 2)
!     Upper part.
m(2, 4) = m(2, 4) / m(2, 2)

!     Lower part.
temp = m(3, 1) * m(1, 3)
temp = temp + m(3, 2) * m(2, 3)
m(3, 3) = m(3, 3) - temp
!     Lower part.
temp = m(4, 1) * m(1, 3)
temp = temp + m(4, 2) * m(2, 3)
m(4, 3) = m(4, 3) - temp
!     Upper part.
temp = m(3, 2) * m(2, 4)
m(3, 4) = ( m(3, 4) - temp ) / m(3, 3)

!     Lower part.
temp = m(4, 2) * m(2, 4)
temp = temp + m(4, 3) * m(3, 4)
m(4, 4) = m(4, 4) - temp

return

!--------
! FORMATS
!--------
#if defined(_CS_LANG_FR)
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : CHIMIE ATMOSPHERIQUE DEMANDEE'                   ,/,&
  '@    ========'                                                ,/,&
  '@    L''utilisateur a choisi de fournir son propre schema'    ,/,&
  '@    chimique. Or aucun schema n''a ete trouve.'              ,/,&
  '@'                                                            ,/,&
  '@  Le calcul ne sera pas execute.'                            ,/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#else
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : ATMOSPHERIC CHEMISTRY'                           ,/,&
  '@    ========'                                                ,/,&
  '@    The user choose to use its own chemical scheme'          ,/,&
  '@    However no scheme has been found'                        ,/,&
  '@'                                                            ,/,&
  '@  The computation will not be run                           ',/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#endif

end subroutine lu_decompose


!===============================================================================

!> lu_solve
!> \brief Resolution of MY=X where M is an LU factorization computed
!>        by lu_decompose
!------------------------------------------------------------------------------

!------------------------------------------------------------------------------
! Arguments
!------------------------------------------------------------------------------
!   mode          name             role
!------------------------------------------------------------------------------
!> \param[in]     ns               matrix row number from the chemical
!>                                 species number
!> \param[in]     m                an LU factorization computed by lu_decompose
!> \param[in,out] x                on entry, the right-hand side of the equation
!                                  on exit, the solution of the equation
!______________________________________________________________________________

subroutine lu_solve (ns, m, x)

use entsor

implicit none

! Arguments

integer ns
double precision m(ns,ns)
double precision x(ns)

! Local variables

double precision temp

if(1.eq.1) then
  write(nfecra,9000)
  call csexit (1)
endif

! Example

!     Forward substitution.

x(1) = x(1) / m(1, 1)

temp = m(2, 1) * x(1)
x(2) = ( x(2) - temp ) / m(2, 2)

temp = m(3, 1) * x(1)
temp = temp + m(3, 2) * x(2)
x(3) = ( x(3) - temp ) / m(3, 3)

temp = m(4, 1) * x(1)
temp = temp + m(4, 2) * x(2)
temp = temp + m(4, 3) * x(3)
x(4) = ( x(4) - temp ) / m(4, 4)


!     Backward substitution.

temp = m(3, 4) * x(4)
x(3) = x(3) - temp

temp = m(2, 3) * x(3)
temp = temp + m(2, 4) * x(4)
x(2) = x(2) - temp

temp = m(1, 3) * x(3)
x(1) = x(1) - temp

return

!--------
! FORMATS
!--------
#if defined(_CS_LANG_FR)
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : CHIMIE ATMOSPHERIQUE DEMANDEE'                   ,/,&
  '@    ========'                                                ,/,&
  '@    L''utilisateur a choisi de fournir son propre schema'    ,/,&
  '@    chimique. Or aucun schema n''a ete trouve.'              ,/,&
  '@'                                                            ,/,&
  '@  Le calcul ne sera pas execute.'                            ,/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#else
  9000 format(                                                           &
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@'                                                            ,/,&
  '@ @@ ERROR : ATMOSPHERIC CHEMISTRY'                           ,/,&
  '@    ========'                                                ,/,&
  '@    The user choose to use its own chemical scheme'          ,/,&
  '@    However no scheme has been found'                        ,/,&
  '@'                                                            ,/,&
  '@  The computation will not be run                           ',/,&
  '@'                                                            ,/,&
  '@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@@',/,&
  '@                                                            ',/)
#endif
end subroutine lu_solve