MFIX  2016-1
parse_mod.f
Go to the documentation of this file.
1 MODULE parse
2 
3  USE compar
4  USE exit, only: mfix_exit
5  USE funits
6  USE param
7  USE param1
8 
9  IMPLICIT NONE
10 
11 ! Strings indicating arithmetic operation and reaction blocks.
12  CHARACTER(LEN=2), PARAMETER :: start_str = '@(' ! start
13  CHARACTER(LEN=1), PARAMETER :: end_str = ')' ! end
14 
15 ! Strings indicating reaction blocks.
16  CHARACTER(LEN=4), PARAMETER :: rxn_blk = 'RXNS' ! start block
17  CHARACTER(LEN=8), PARAMETER :: des_rxn_blk = 'DES_RXNS' ! start block
18  CHARACTER(LEN=3), PARAMETER :: end_blk = 'END' ! end block
19 
20  LOGICAL reading_rxn
21  LOGICAL reading_rate
22 
23  LOGICAL des_rxn
24  LOGICAL tfm_rxn
25 
26 ! Logical indicating that the start of a reaction construct has
27 ! been identified.
28  LOGICAL in_construct
29 
30 ! Logical indicating that the chemical equation spans additional lines.
31  LOGICAL more_chemeq
32 
33 ! Reaction names
34  CHARACTER(len=32), DIMENSION(:), ALLOCATABLE :: rxn_name
35 ! chemical Equations
36  CHARACTER(len=512), DIMENSION(:), ALLOCATABLE :: rxn_chem_eq
37 ! User defined heat of reaction
38  DOUBLE PRECISION, DIMENSION(:), ALLOCATABLE :: usrdh
39 ! User defined heat of reaction partitions.
40  DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE :: usrfdh
41 
42 
43 ! Logical indicating that the start of a reaction construct has
44 ! been identified.
46 
47 ! Reaction names
48  CHARACTER(len=32), DIMENSION(:), ALLOCATABLE :: des_rxn_name
49 ! chemical Equations
50  CHARACTER(len=512), DIMENSION(:), ALLOCATABLE :: des_rxn_chem_eq
51 ! User defined heat of reaction
52  DOUBLE PRECISION, DIMENSION(:), ALLOCATABLE :: des_usrdh
53 ! User defined heat of reaction partitions.
54  DOUBLE PRECISION, DIMENSION(:,:), ALLOCATABLE :: des_usrfdh
55 
56  CONTAINS
57 
58 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
59 ! Function name: setReaction !
60 ! !
61 ! Purpose: !
62 ! !
63 ! Variables referenced: None !
64 ! !
65 ! Variables modified: None !
66 ! !
67 ! Local variables: None !
68 ! !
69 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
70  SUBROUTINE setreaction(RxN, lNg, lSAg, lM, lNs, lSAs, lDH, lfDH)
71 
72  use rxn_com
73  use toleranc
74 
75  IMPLICIT NONE
76 
77 ! Pass Arguments:
78 !---------------------------------------------------------------------//
79 ! Data structure for storing reaction data.
80  TYPE(reaction_block), POINTER, INTENT(INOUT) :: RxN
81 ! Number of gas speices
82  INTEGER, INTENT(IN) :: lNg
83 ! Gas phase species aliases
84  CHARACTER(len=32), DIMENSION(DIM_N_g), INTENT(IN) :: lSAg
85 ! Number of solids phases
86  INTEGER, INTENT(IN) :: lM
87 ! Number of species in each solids phase.
88  INTEGER, DIMENSION(DIM_M), INTENT(IN) :: lNs
89 ! Solids phase speices aliases.
90  CHARACTER(len=32), DIMENSION(DIM_M, DIM_N_s), INTENT(IN) :: lSAs
91 ! User defined heat of reaction.
92  DOUBLE PRECISION, INTENT(IN) :: lDH
93 ! User defined heat of reaction partition.
94  DOUBLE PRECISION, DIMENSION(0:DIM_M), INTENT(IN) :: lfDH
95 
96 ! Local Variables:
97 !---------------------------------------------------------------------//
98 ! Alias, phase, species, stoich coeff :: Reactants Products
99  CHARACTER(LEN=32), DIMENSION(50) :: rAlias , pAlias
100  DOUBLE PRECISION, DIMENSION(50) :: rCoeff , pCoeff
101 
102 ! Number of products and reactants
103  INTEGER rNo, pNo
104 ! Positions in ChemEq delineating the reactants and products.
105  INTEGER rEnd, pStart
106 ! Loop counters
107  INTEGER L, LL, M, lN
108 
109 ! Sum of user specified heat of reaction partitions. If fracDH is set
110 ! by the user, they must sum to one over all phases.
111  DOUBLE PRECISION sumFDH
112 ! Local storage for chemical equations, left-adjusted and trimmed
113  CHARACTER(LEN=512) lChemEq
114 ! Local storage for reaction name, left-adjusted and trimmed
115  CHARACTER(LEN=32) lName
116 ! Logical indicating the reaction is skipped.
117  LOGICAL Skip
118 
119  LOGICAL pMap(0:lm)
120 
121  INTEGER nSpecies, nPhases
122 
123  LOGICAL blankAlias(0:(dim_n_g + lm*dim_n_s))
124 
125 ! Initialize local reaction name and chemical equation variables.
126  lname = trim(adjustl(rxn%Name))
127  lchemeq = trim(adjustl(rxn%ChemEq))
128 
129  rxn%Classification = "Undefined"
130  rxn%Calc_DH = .true.
131  rxn%nSpecies = 0
132  rxn%nPhases = 0
133 
134 ! Verify that the reactants are separated by --> or = signs. If the
135 ! chemical equation is NONE, the reaction is skipped.
136  CALL checksplit(lname, lchemeq, rend, pstart, skip)
137  IF(skip) THEN
138  rxn%nSpecies = 0
139  RETURN
140  ENDIF
141 ! Set the flag to calculate heat of reaction.
142  rxn%Calc_DH = .true.
143  IF(ldh /= undefined) rxn%Calc_DH = .false.
144 
145 ! Pull off the reactants from the chemical equations.
146  CALL splitentries(lname, lchemeq, 1, rend, rno, ralias, rcoeff)
147 ! Pull off the products from the chemical equations.
148  CALL splitentries(lname, lchemeq, pstart, len_trim(lchemeq), &
149  pno, palias, pcoeff)
150 
151  nspecies = rno + pno
152  rxn%nSpecies = nspecies
153  Allocate( rxn%Species( nspecies ))
154 
155  CALL checkblankaliases(lng, lsag, lm, lns, lsas, blankalias)
156 
157 ! Check that species in the chemical equation match a species alias
158 ! in one of the phases.
159  CALL mapaliases(lname, lchemeq, lng, lsag, lm, lns, lsas, rno, &
160  ralias, rcoeff, -one, 0, blankalias, rxn)
161 
162 ! Check that species in the chemical equation match a species alias
163 ! in one of the phases.
164  CALL mapaliases(lname, lchemeq, lng, lsag, lm, lns, lsas, pno, &
165  palias, pcoeff, one, rno, blankalias, rxn)
166 
167 ! All the relevant data has been collected at this point. Build the
168 ! reaction block data structure.
169  l = max(1,lm)
170  ll = (l * (l-1)/2)
171  Allocate( rxn%rPhase( ll+l ))
172 
173 
174 ! Initialize local map and global values
175  pmap(:) = .false.
176  nphases = 0
177  DO ln = 1, nspecies
178  m = rxn%Species(ln)%pMap
179 
180  rxn%Species(ln)%mXfr = m
181  rxn%Species(ln)%xXfr = zero
182 
183  IF(.NOT.pmap(m)) THEN
184  pmap(m) = .true.
185  nphases = nphases + 1
186  ENDIF
187  ENDDO
188  rxn%nPhases = nphases
189 
190 ! Initialize sum of heat of reaction partitions.
191  sumfdh = zero
192 ! The user specified the heat of reaction.
193  IF(.NOT.rxn%Calc_DH) THEN
194 ! Allocate and initialize the heat of reaction storage array.
195  Allocate( rxn%HoR( 0:lm ))
196  rxn%HoR(:) = zero
197  DO m=0,lm
198 ! The phase is referenced by the reaction and heat of reaction is
199 ! allocated (in part or fully) this this phase.
200  IF(pmap(m) .AND. lfdh(m) .NE. undefined) THEN
201 ! Store the heat of reaction.
202  rxn%HoR(m) = lfdh(m) * ldh
203  sumfdh = sumfdh + lfdh(m)
204 ! The phase is not referenced by the reaction, but the heat of reaction
205 ! is allocated (in part or fully) to this phase. Flag error and exit.
206  ELSEIF(.NOT.pmap(m) .AND. lfdh(m) .NE. undefined) THEN
207  IF(mype == pe_io) THEN
208  write(*,1000) trim(lname)
209  write(*,1001)
210  write(*,1010)
211 ! Write message to log file.
212  write(unit_log,1000) trim(lname)
213  write(unit_log,1001)
214  write(unit_log,1010)
215  ENDIF
216 ! Exit MFIX
217  CALL mfix_exit(mype)
218  ENDIF
219  ENDDO
220 ! Logical check: No partition was assigned to an undefined phase.
221  DO m=lm+1,dim_m
222  IF(.NOT.rxn%Calc_DH .AND. lfdh(m) .NE. undefined) THEN
223  IF(mype == pe_io) THEN
224  write(*,1000) trim(lname)
225  write(*,1001)
226  write(*,1010)
227 ! Write message to log file.
228  write(unit_log,1000) trim(lname)
229  write(unit_log,1001)
230  write(unit_log,1010)
231 ! Exit MFIX
232  ENDIF
233  CALL mfix_exit(mype)
234  ENDIF
235  ENDDO
236  ENDIF
237 
238 ! Verify that the heat of reaction partitions sum to one.
239  IF(.NOT.rxn%Calc_DH .AND. .NOT. compare(sumfdh, one)) THEN
240  IF(mype == pe_io) THEN
241  write(*,1002) trim(lname)
242  write(*,1010)
243 ! Write message to log file.
244  write(unit_log,1002) trim(lname)
245  write(unit_log,1010)
246 ! Exit MFIX
247  CALL mfix_exit(mype)
248  ENDIF
249  ENDIF
250 
251  RETURN
252 
253  1000 FORMAT(/1x,70('*')/' From: From: setReaction:',/ &
254  ' Message: Heat of reaction is proportioned to a phase not', &
255  ' referenced',/' by the chemical equation for reaction ',a,'.')
256 
257  1001 FORMAT(/' If this is a catalytic reaction, reference one of the',&
258  ' species of the',/' catalyst phase within the chemical', &
259  ' equation with a stoichiometric',/' coefficient of zero.'/)
260 
261  1002 FORMAT(/1x,70('*')/' From: From: setReaction:',/ &
262  ' Message: The heat of reaction partitions (fracDH) to all', &
263  ' phases do',/' not sum to one for reaction ',a,'.')
264 
265  1010 FORMAT(' Please refer to the Readme file for chemical equation', &
266  ' input formats',/' and correct the data file.',/1x,70('*')/)
267 
268 
269  END SUBROUTINE setreaction
270 
271 
272 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
273 ! Function name: checkSplit !
274 ! !
275 ! Purpose: Determine the location of reactatns and products within !
276 ! the chemical equation. If the entry is NONE, flag that the reaction !
277 ! is to be skipped for further processing. !
278 ! !
279 ! Variables referenced: None !
280 ! !
281 ! Variables modified: None !
282 ! !
283 ! Local variables: None !
284 ! !
285 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
286  SUBROUTINE checksplit(lName, lChemEq, lrEnd, lpStart, lSkip)
288  IMPLICIT NONE
289 
290 ! Pass Arguments:
291 !---------------------------------------------------------------------//
292 ! Chemical reaction name.
293  CHARACTER(len=*), INTENT(IN) :: lName
294 ! Chemical Equation from deck file.
295  CHARACTER(len=*), INTENT(IN) :: lChemEq
296 ! Position specifying the end of the reactants in lChemEq
297  INTEGER, INTENT(OUT) :: lrEnd
298 ! Position specifying the start of the products in lChemEq
299  INTEGER, INTENT(OUT) :: lpStart
300 ! If the chemical equation is NONE, the the skip flag is set to avoid
301 ! further processing.
302  LOGICAL, INTENT(OUT) :: lSkip
303 
304 ! Local Variables:
305 !---------------------------------------------------------------------//
306 ! Position of the head (>) and tail (-) of an arror (-->)
307  INTEGER hArr, tArr
308 ! Postion of the head and tail of equal signs (=, ==, ===, ...)
309  INTEGER hEqs, tEqs
310 ! Position of the head/tail of a reverse arrow (<--)
311  INTEGER hRArr, tRArr
312 ! A flag generated to point out the location of the entry error.
313  CHARACTER(LEN=512) FLAG
314  flag = ''
315 
316 ! If the chemical equation is set to 'none', then the reaction is
317 ! skipped for the simulation.
318  lskip = .false.
319  IF(index(lchemeq,'NONE') > 0) THEN
320  lskip = .true.
321  lrend = undefined_i
322  lpstart = undefined_i
323  RETURN
324  ENDIF
325 
326 ! Search for > (part of -->) and search for < (part of <--)
327  tarr = index(lchemeq,'-', back=.false.)
328  harr = index(lchemeq,">", back=.true.)
329 ! Search for the first and last instances of equal signs.
330  teqs = index(lchemeq,"=", back=.false.)
331  heqs = index(lchemeq,"=", back=.true.)
332 ! Search for < (as part of <-- or <-->). Illegal chemical equation.
333  hrarr = index(lchemeq,"<", back=.false.)
334  trarr = index(lchemeq,"-", back=.true.)
335 
336 ! An illegal arrow was found! Flag error and exit.
337  IF(hrarr > 0) THEN
338  IF(mype == pe_io) THEN
339 ! Construct the error flag.
340  IF(harr > 0) THEN
341  flag = setflag(20, hrarr, harr)
342  ELSEIF(trarr > 0) THEN
343  flag = setflag(20, hrarr, trarr)
344  ELSE
345  flag = setflag(20, hrarr)
346  ENDIF
347 ! Write the message to the screen.
348  write(*,1000) trim(lname)
349  write(*,1002)'Illegal'
350  write(*,1010) trim(lchemeq), trim(flag)
351  write(*,1001)
352 ! Write message to log file.
353  write(unit_log,1000) trim(lname)
354  write(unit_log,1002)'Illegal'
355  write(unit_log,1010) trim(lchemeq), trim(flag)
356  write(unit_log,1001)
357 ! Exit MFIX
358  CALL mfix_exit(mype)
359  ENDIF
360  ENDIF
361 ! If there are more than one operator, flag error and exit.
362  IF(harr /= 0 .AND. heqs /= 0) THEN
363  IF(mype == pe_io) THEN
364 ! Construct the error flag.
365  flag = setflag(20, harr, heqs)
366 ! Write the message to the screen.
367  write(*,1000) trim(lname)
368  write(*,1002)'Too many'
369  write(*,1010) trim(lchemeq), trim(flag)
370  write(*,1001)
371 ! Write message to log file.
372  write(unit_log,1000) trim(lname)
373  write(unit_log,1002)'Too many'
374  write(unit_log,1010) trim(lchemeq), trim(flag)
375  write(unit_log,1001)
376 ! Exit MFIX
377  CALL mfix_exit(mype)
378  ENDIF
379 ! If there is no operator (--> or =), flag error and exit.
380  ELSEIF(harr == 0 .AND. heqs == 0) THEN
381  IF(mype == pe_io) THEN
382  write(*,1000) trim(lname)
383  write(*,1002) 'No'
384  write(*,1011) trim(lchemeq)
385  write(*,1001)
386 ! Write message to log file.
387  write(unit_log,1000) trim(lname)
388  write(unit_log,1002) 'No'
389  write(unit_log,1011) trim(lchemeq)
390  write(unit_log,1001)
391 ! Exit MFIX
392  CALL mfix_exit(mype)
393  ENDIF
394 ! The head of an arrow was found.
395  ELSEIF(harr /= 0) THEN
396 ! Verify that a tail was found.
397  IF(tarr == 0) THEN
398 ! Construct the error flag.
399  flag = setflag(20, harr)
400  IF(mype == pe_io) THEN
401  write(*,1000) trim(lname)
402  write(*,1003) 'Missing the tail; -->'
403  write(*,1010) trim(lchemeq), trim(flag)
404  write(*,1001)
405 ! Write message to log file.
406  write(unit_log,1000) trim(lname)
407  write(unit_log,1003) 'Missing the tail; -->'
408  write(unit_log,1010) trim(lchemeq), trim(flag)
409  write(unit_log,1001)
410 ! Exit MFIX
411  CALL mfix_exit(mype)
412  ENDIF
413  ELSEIF(tarr > harr) THEN
414  IF(mype == pe_io) THEN
415  flag = setflag(20, harr, index(lchemeq,'-',back=.true.))
416  write(*,1000) trim(lname)
417  write(*,1003) 'Arror head preceeds the tail; -->'
418  write(*,1010) trim(lchemeq), trim(flag)
419  write(*,1001)
420 ! Write message to log file.
421  write(unit_log,1000) trim(lname)
422  write(unit_log,1003) 'Arror head preceeds the tail; -->'
423  write(unit_log,1010) trim(lchemeq), trim(flag)
424  write(unit_log,1001)
425 ! Exit MFIX
426  CALL mfix_exit(mype)
427  ENDIF
428  ELSE
429 ! An arror was used to seperate reactants and products. Send back the
430 ! ending index of reactants and the starting index for products.
431  lrend = tarr - 1
432  lpstart = harr + 1
433  ENDIF
434 ! Equals sign(s) were used to specify the reaction. Send back the ending
435 ! index of reactants and the starting index for products.
436  ELSEIF(heqs /= 0) THEN
437  lrend = teqs - 1
438  lpstart = heqs + 1
439 ! Fatal Error. One of the above checks should have caught any problems
440 ! and sent out an error message.
441  ELSE
442  IF(mype == pe_io) THEN
443  write(*,1000) trim(lname)
444  write(*,1004)
445  write(*,1011) trim(lchemeq)
446  write(*,1001)
447 ! Write message to log file.
448  write(unit_log,1000) trim(lname)
449  write(unit_log,1004)
450  write(unit_log,1010) trim(lchemeq), trim(flag)
451  write(unit_log,1001)
452 ! Exit MFIX
453  CALL mfix_exit(mype)
454  ENDIF
455  ENDIF
456 
457  RETURN
458 
459  1000 FORMAT(/1x,70('*')/' From: From: setReaction --> checkSplit',/ &
460  ' Message: Error in determining the reactants and products', &
461  ' in the',/' chemical equation for reaction ',a,'.')
462 
463  1001 FORMAT(' Please refer to the Readme file for chemical equation', &
464  ' input formats',/' and correct the data file.',/1x,70('*')/)
465 
466  1002 FORMAT(/1x,a,' operators were found!')
467 
468  1003 FORMAT(' Incorrect operator format! ',a)
469 
470  1004 FORMAT(' FATAL ERROR: All logical checks failed.')
471 
472  1010 FORMAT(/' Chemical Equation: ',a,/1x, a/)
473 
474  1011 FORMAT(/' Chemical Equation: ',a,/)
475 
476 
477  END SUBROUTINE checksplit
478 
479 
480 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
481 ! Function name: splitEntries() !
482 ! !
483 ! Purpose: Takes a string of either reactants or products and splits !
484 ! the string into individual species/stoichiometric entries. !
485 ! !
486 ! A call to splitAliasAndCoeff is made to further split the entries !
487 ! into species aliases and matching stoichiometric coefficients. !
488 ! !
489 ! Variables referenced: None !
490 ! !
491 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
492  SUBROUTINE splitentries(lName, lChemEq, lStart, lEnd, lNo, &
493  lalias, lcoeff)
495  IMPLICIT NONE
496 
497 ! Pass Arguments:
498 !---------------------------------------------------------------------//
499 ! Chemical reaction name.
500  CHARACTER(len=*), INTENT(IN) :: lName
501 ! Chemical Equation.
502  CHARACTER(len=*), INTENT(IN) :: lChemEq
503 ! Starting position for substring analysis.
504  INTEGER, INTENT(IN) :: lStart
505 ! Ending position for substring analysis.
506  INTEGER, INTENT(IN) :: lEnd
507 ! The number of individual species found in lSpecies.
508  INTEGER, INTENT(OUT) :: lNo
509 ! Species Aliases from the chemical equation.
510  CHARACTER(LEN=32), DIMENSION(50), INTENT(OUT) :: lAlias
511 ! Stoichiometric coefficient pulled from the chemical equation.
512  DOUBLE PRECISION, DIMENSION(50), INTENT(OUT) :: lCoeff
513 
514 ! Local Variables:
515 !---------------------------------------------------------------------//
516 ! Flag indicating that there are more entries to process.
517  LOGICAL MORE
518 ! Starting position for left-to-right search.
519  INTEGER lPOS
520 ! Position of plus sign charachter found in search.
521  INTEGER rPOS
522 
523 ! Initialize storage variables.
524  lno = 0
525  lalias(:) = ''
526  lcoeff(:) = undefined
527 
528 ! Initialiase local variables.
529  lpos = lstart
530  more = .true.
531 ! Loop through the string, splitting entries separated by a plus sign.
532  DO WHILE(more)
533 ! Increment the species counter.
534  lno = lno + 1
535 ! Locate the plus sign. (Left to right)
536  rpos = (lpos-1) + index(lchemeq(lpos:lend),"+", back=.false.)
537 ! A plus sign was found.
538  IF(rpos .GT. lpos) THEN
539 ! Extract the entry and split it into the species alias and
540 ! stoichiometric coefficient.
541  CALL splitaliasandcoeff(lname, lchemeq, lpos, rpos-1, &
542  lalias(lno), lcoeff(lno))
543 ! Indicate that there are more entries to process.
544  more = .true.
545 ! No plus sign was found. This is the last entry.
546  ELSE
547 ! Extract the entry and split it into the species alias and
548 ! stoichiometric coefficient.
549  CALL splitaliasandcoeff(lname, lchemeq, lpos, lend, &
550  lalias(lno), lcoeff(lno))
551 ! Indicate that there are no more entries to process.
552  more = .false.
553  ENDIF
554 ! Move past the found plus sign for next search.
555  lpos = rpos + 1
556  ENDDO
557 
558  RETURN
559  END SUBROUTINE splitentries
560 
561 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
562 ! Function name: splitAliasAndCoeff() !
563 ! !
564 ! Purpose: Take a string containing a species alias and stoichio- !
565 ! metric coefficient and splits them into their respective parts. !
566 ! !
567 ! If no numerical coefficient is found, it is set to one. !
568 ! !
569 ! If present, asterisks (*) are assumed to seperate numerical !
570 ! coefficients and the species alias. If more than one asterisk is !
571 ! found, and error is reported and MFIX_EXIT is called. !
572 ! !
573 ! Variables referenced: None !
574 ! !
575 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
576  SUBROUTINE splitaliasandcoeff(lName, lChemEq, lStart, lEnd, &
577  lalias, lcoeff)
579  IMPLICIT NONE
580 
581 ! Pass Arguments:
582 !---------------------------------------------------------------------//
583 ! Chemical reaction name.
584  CHARACTER(len=*), INTENT(IN) :: lName
585 ! Chemical Equation.
586  CHARACTER(len=*), INTENT(IN) :: lChemEq
587 ! Starting position for substring analysis.
588  INTEGER, INTENT(IN) :: lStart
589 ! Ending position for substring analysis.
590  INTEGER, INTENT(IN) :: lEnd
591 ! Species Aliases from the chemical equation.
592  CHARACTER(LEN=32), INTENT(OUT) :: lAlias
593 ! Stoichiometric coefficient pulled from the chemical equation.
594  DOUBLE PRECISION, INTENT(OUT) :: lCoeff
595 
596 ! Local Variables:
597 !---------------------------------------------------------------------//
598 ! Flag
599  LOGICAL MATCH
600  INTEGER nPOS
601 
602  INTEGER L, N, IOS, aPOS, a2POS
603 
604  CHARACTER(LEN=12), PARAMETER :: Numbers = '.0123456789'
605 
606 ! A flag generated to point out the location of the entry error.
607  CHARACTER(LEN=512) FLAG
608  flag = ''
609 
610 ! Locate the first asterisk (if any). Search left-to-right.
611  apos = index(lchemeq(lstart:lend),"*", back=.false.)
612 ! An asterisk was found.
613  IF(apos .GT. zero) THEN
614 ! Make sure that there isn't another asterisk further down the string.
615  a2pos = index(lchemeq(lstart:lend),"*", back=.true.)
616  IF(apos /= a2pos) THEN
617  IF(mype == pe_io) THEN
618 ! Construct the error flag.
619  flag = setflag(20, lstart+apos, lstart+a2pos)
620 ! Write the message to the screen.
621  write(*,1000) trim(lname)
622  write(*,1002)'Too many'
623  write(*,1010) trim(lchemeq), trim(flag)
624  write(*,1001)
625 ! Write message to log file.
626  write(unit_log,1000) trim(lname)
627  write(unit_log,1002)'Too many'
628  write(unit_log,1010) trim(lchemeq), trim(flag)
629  write(unit_log,1001)
630 ! Exit MFIX
631  CALL mfix_exit(mype)
632  ENDIF
633  ELSE
634 ! Store left-of-asterisk as the coefficient. If an error occurs in
635 ! converting the string to double precision, flag the problem and
636 ! call MFIX_EXIT.
637  READ(lchemeq(lstart:(lstart+apos-2)),*,iostat=ios) lcoeff
638  IF(ios .NE. 0 .AND. mype == pe_io) THEN
639 ! Construct the error flag.
640  flag = setflag(20, lstart + int(apos/2))
641 ! Write the message to the screen.
642  write(*,1000) trim(lname)
643  write(*,1010) trim(lchemeq), trim(flag)
644  write(*,1001)
645 ! Write message to log file.
646  write(unit_log,1000) trim(lname)
647  write(unit_log,1010) trim(lchemeq), trim(flag)
648  write(unit_log,1001)
649 ! Exit MFIX
650  CALL mfix_exit(mype)
651  ENDIF
652 ! Store right-of-asterisk as the species alias.
653  WRITE(lalias,"(A)") &
654  trim(adjustl(lchemeq((lstart+apos):lend)))
655  ENDIF
656 ! If no asterisk was found, search for numbers and spaces.
657  ELSE
658 ! Initialize the position of last consecutive number.
659  npos = 0
660 ! In a left-to-right search, check if the characters in the entry are
661 ! numbers or punctuation.
662  DO l=lstart,lend
663  match = .false.
664  DO n=1,12
665  IF(lchemeq(l:l) /= numbers(n:n)) cycle
666 ! Note the position of the number.
667  npos = l
668 ! Flag that a match was made.
669  match = .true.
670  ENDDO
671 ! If no match, assume the end of the coefficient was found.
672  IF(.NOT.match) EXIT
673  ENDDO
674 ! If no numbers or punctuation was found, assumed the stoichiometric
675 ! coefficient is one.
676  IF(trim(lchemeq(lstart:npos)) =='') THEN
677  lcoeff = 1.0d0
678  ELSE
679 ! If leading numbers were found, store as the stoich-coeff.
680  READ(lchemeq(lstart:npos),*,iostat=ios) lcoeff
681 ! Report any problems in converting the string to double precision.
682  IF(ios .NE. 0 .AND. mype == pe_io) THEN
683 ! Construct the error flag.
684  flag = setflag(20, &
685  lstart+int(len_trim(lchemeq(lstart:npos))/2))
686 ! Write the message to the screen.
687  write(*,1000) trim(lname)
688  write(*,1010) trim(lchemeq), trim(flag)
689  write(*,1001)
690 ! Write message to log file.
691  write(unit_log,1000) trim(lname)
692  write(unit_log,1010) trim(lchemeq), trim(flag)
693  write(unit_log,1001)
694 ! Exit MFIX
695  CALL mfix_exit(mype)
696  ENDIF
697  ENDIF
698 ! Store right-of-coefficient as the species alias.
699  READ(lchemeq(npos+1:lend),*,iostat=ios) lalias
700  ENDIF
701 ! Quick check to make sure that the species alias is not empty.
702  IF(len_trim(lalias) == 0 .AND. mype == pe_io) THEN
703 ! Construct the error flag.
704  flag = setflag(20, lstart + int(lend/2))
705 ! Write the message to the screen.
706  write(*,1003) trim(lname)
707  write(*,1010) trim(lchemeq), trim(flag)
708  write(*,1001)
709 ! Write message to log file.
710  write(unit_log,1003) trim(lname)
711  write(unit_log,1010) trim(lchemeq), trim(flag)
712  write(unit_log,1001)
713 ! Exit MFIX
714  CALL mfix_exit(mype)
715  ENDIF
716 
717  RETURN
718 
719  1000 FORMAT(/1x,70('*')/' From: From: setReaction -->', &
720  ' splitAliasAndCoeff',/' Message: Error determining the', &
721  ' stoichiometric coefficient in the',/' chemical equation', &
722  ' for reaction ',a,'.')
723 
724 
725  1001 FORMAT(' Please refer to the Readme file for chemical equation', &
726  ' input formats',/' and correct the data file.',/1x,70('*')/)
727 
728  1002 FORMAT(/1x,a,' operators were found!')
729 
730  1003 FORMAT(/1x,70('*')/' From: From: setReaction -->', &
731  ' splitAliasAndCoeff',/' Message: Error determining the', &
732  ' speices in the chemical equation for',/' reaction ',a,'.'/)
733 
734  1010 FORMAT(/' Chemical Equation: ',a,/1x, a/)
735 
736  1011 FORMAT(/' Chemical Equation: ',a,/)
737 
738  END SUBROUTINE splitaliasandcoeff
739 
740 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
741 ! Function name: checkBlankAliases() !
742 ! !
743 ! Purpose: Take a string containing a species alias and stoichio- !
744 ! metric coefficient and splits them into their respective parts. !
745 ! !
746 ! Variables referenced: None !
747 ! !
748 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
749  SUBROUTINE checkblankaliases(lNg, lSAg, lM, lNs, lSAs, lBA)
751  IMPLICIT NONE
752 
753 ! Number of gas speices
754  INTEGER, INTENT(IN) :: lNg
755 ! Gas phase species aliases
756  CHARACTER(len=32), DIMENSION(DIM_N_g), INTENT(IN) :: lSAg
757 ! Number of solids phases
758  INTEGER, INTENT(IN) :: lM
759 ! Number of species in each solids phase.
760  INTEGER, DIMENSION(DIM_M), INTENT(IN) :: lNs
761 ! Solids phase speices aliases.
762  CHARACTER(len=32), DIMENSION(DIM_M, DIM_N_s), INTENT(IN) :: lSAs
763 
764  LOGICAL, INTENT(OUT) :: lBA(0:(dim_n_g + lm*dim_n_s))
765 
766  INTEGER M, N
767 
768 ! Loop counter for continuum and discrete species
769  INTEGER Nsp
770 
771 ! Initialize counters
772  nsp = 0
773 
774  lba(0) = .false.
775  DO n = 1, lng
776  nsp = nsp + 1
777  lba(nsp) = .false.
778  IF(len_trim(lsag(n)) == 0) THEN
779  lba(nsp) = .true.
780  lba(0) = .true.
781  ENDIF
782  ENDDO
783 
784 ! Compare aliaes between solids phases
785  DO m = 1, lm
786  DO n = 1, lns(m)
787  nsp = nsp + 1
788  lba(nsp) = .false.
789  IF(len_trim(lsas(m,n)) == 0) THEN
790  lba(nsp) = .true.
791  lba(0) = .true.
792  ENDIF
793  ENDDO
794  ENDDO
795 
796  RETURN
797  END SUBROUTINE checkblankaliases
798 
799 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
800 ! Function name: mapAliases() !
801 ! !
802 ! Purpose: Take a string containing a species alias and stoichio- !
803 ! metric coefficient and splits them into their respective parts. !
804 ! !
805 ! If no numerical coefficient is found, it is set to one. !
806 ! !
807 ! If present, asterisks (*) are assumed to seperate numerical !
808 ! coefficients and the species alias. If more than one asterisk is !
809 ! found, and error is reported and MFIX_EXIT is called. !
810 ! !
811 ! Variables referenced: None !
812 ! !
813 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
814  SUBROUTINE mapaliases(lName, lChemEq, lNg, lSAg, lM, lNs, lSAs, &
815  lno, lalias, lcoeff, lsgn, lstart, lba, lrxn)
817  use rxn_com
818 
819  IMPLICIT NONE
820 
821 ! Pass Arguments:
822 !---------------------------------------------------------------------//
823 ! Chemical reaction name.
824  CHARACTER(len=*), INTENT(IN) :: lName
825 ! Chemical Equation.
826  CHARACTER(len=*), INTENT(IN) :: lChemEq
827 ! Number of gas speices
828  INTEGER, INTENT(IN) :: lNg
829 ! Gas phase species aliases
830  CHARACTER(len=32), DIMENSION(DIM_N_g), INTENT(IN) :: lSAg
831 ! Number of solids phases
832  INTEGER, INTENT(IN) :: lM
833 ! Number of species in each solids phase.
834  INTEGER, DIMENSION(DIM_M), INTENT(IN) :: lNs
835 ! Solids phase speices aliases.
836  CHARACTER(len=32), DIMENSION(DIM_M, DIM_N_s), INTENT(IN) :: lSAs
837 ! Number of products (or reactants)
838  INTEGER, INTENT(IN) :: lNo
839 ! Species Alaises pulled from the chemical equation.
840  CHARACTER(LEN=32), DIMENSION(50), INTENT(IN) :: lAlias
841 
842  DOUBLE PRECISION, DIMENSION(50), INTENT(IN) :: lCoeff
843 
844  DOUBLE PRECISION, INTENT(IN) :: lSgn
845 
846  INTEGER, INTENT(IN) :: lStart
847 
848  LOGICAL, INTENT(IN) :: lBA(0:(dim_n_g + lm*dim_n_s))
849 
850 ! Data structure for storing reaction data.
851  TYPE(reaction_block), POINTER, INTENT(INOUT) :: lRxN
852 
853 
854 ! Local Variables:
855 !---------------------------------------------------------------------//
856 ! Loop counters.
857  INTEGER L, M, NN
858 ! A flag generated to point out the location of the entry error.
859  CHARACTER(LEN=512) FLAG
860 ! Location in string to locate error.
861  INTEGER lPOS, rPOS
862 
863 ! Initialize the error flag.
864  flag = ''
865 
866 ! Loop over the number of (reactants/products)
867  aloop : DO l=1, lno
868 ! Compare entry with gas phase species.
869  DO nn = 1, lng
870  IF( checkmatch(lsag(nn), lalias(l))) THEN
871  lrxn%Species(lstart + l)%pMap = 0
872  lrxn%Species(lstart + l)%sMap = nn
873  lrxn%Species(lstart + l)%Coeff = lsgn * lcoeff(l)
874  cycle aloop
875  ENDIF
876  ENDDO
877 ! Compare entry with solids phase species.
878  DO m = 1, lm
879  DO nn = 1, lns(m)
880  IF(checkmatch(lsas(m,nn),lalias(l))) THEN
881  lrxn%Species(lstart + l)%pMap = m
882  lrxn%Species(lstart + l)%sMap = nn
883  lrxn%Species(lstart + l)%Coeff = lsgn * lcoeff(l)
884  cycle aloop
885  ENDIF
886  ENDDO
887  ENDDO
888 ! No matching species was located. Flag an error and exit.
889  IF(mype == pe_io) THEN
890 ! Construct the error flag.
891  lpos = index(lchemeq,trim(lalias(l)), back=.false.)
892  rpos = index(lchemeq,trim(lalias(l)), back=.true.)
893  flag = setflag(20, 1 + int((lpos + rpos)/2))
894 ! Write the message to the screen.
895  write(*,1000) trim(lalias(l)), trim(lname)
896  write(*,1010) trim(lchemeq), trim(flag)
897  IF(lba(0)) CALL writeba()
898  write(*,1001)
899 ! Write message to log file.
900  write(unit_log,1000) trim(lalias(l)), trim(lname)
901  write(unit_log,1010) trim(lchemeq), trim(flag)
902  IF(lba(0)) CALL writeba(unit_log)
903  write(unit_log,1001)
904  ENDIF
905 ! Exit MFIX
906  CALL mfix_exit(mype)
907 
908  ENDDO aloop
909 
910  RETURN
911 
912  1000 FORMAT(/1x,70('*')/' From: From: setReaction --> mapAliases',/ &
913  ' Message: Unable to match species ',a,' in the chemical', &
914  ' equation for ',/' reaction ',a,'.')
915 
916  1001 FORMAT(/' Please refer to the Readme file for chemical equation',&
917  ' input formats',/' and correct the data file.',/1x,70('*')/)
918 
919  1010 FORMAT(/' Chemical Equation: ',a,/1x, a/)
920 
921  1011 FORMAT(/' Chemical Equation: ',a,/)
922 
923 
924  contains
925 
926 !......................................................................!
927 ! Function name: checkMatch !
928 ! !
929 ! Purpose: Takes two species aliases as arguments, converts them to !
930 ! uppercase and checks if they match. !
931 ! !
932 ! Variables referenced: None !
933 ! !
934 !......................................................................!
935  LOGICAL FUNCTION checkmatch(lSA, ceSA)
937  IMPLICIT NONE
938 
939 ! Pass Arguments:
940 !---------------------------------------------------------------------//
941  CHARACTER(LEN=32), INTENT(IN) :: lSA, ceSA
942 
943 ! Local Variables:
944 !---------------------------------------------------------------------//
945  CHARACTER(LEN=32) tlSA
946 
947 ! Copy species alias.
948  tlsa = lsa
949 ! Remove case sensitivity.
950  CALL make_upper_case (tlsa,32)
951 ! Compare the two strings.
952  checkmatch = .false.
953  IF(trim(tlsa) == trim(cesa)) checkmatch = .true.
954  RETURN
955  END FUNCTION checkmatch
956 
957 !......................................................................!
958 ! Function name: updateMap !
959 ! !
960 ! Purpose: Flags that the passed phase is part of the chemical !
961 ! reaction. If the phase was not already noted, the number of phases !
962 ! in the reaction is increased and the flag set true. Additionally, !
963 ! The number of species (either product or reactant) for the phase !
964 ! is incremented. !
965 ! !
966 ! Variables referenced: None !
967 ! !
968 !......................................................................!
969  SUBROUTINE updatemap(lnP, lpMap, llNoP)
971  IMPLICIT NONE
972 
973 ! Pass Arguments:
974 !---------------------------------------------------------------------//
975  INTEGER, INTENT(INOUT) :: lnP
976 ! Map of phases for the reaction.
977  LOGICAL, INTENT(INOUT) :: lpMap
978  INTEGER, INTENT(INOUT) :: llNoP
979 
980 ! Local Variables:
981 !---------------------------------------------------------------------//
982 ! None
983 
984 ! Increment the number of reactants/products this phase has involved in
985 ! the current reaction.
986  llnop = llnop + 1
987 ! If the phase was already identified, return.
988  IF(lpmap) RETURN
989 ! If this is the first time the phase is identififed, set the flag to
990 ! true and increment the total number of phases in the reaction.
991  lnp = lnp + 1
992  lpmap = .true.
993 
994  RETURN
995  END SUBROUTINE updatemap
996 
997 
998 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
999 ! Function name: writeBA() !
1000 ! !
1001 ! Purpose: Print out which species were not given aliases. !
1002 ! !
1003 ! Variables referenced: None !
1004 ! !
1005 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
1006  SUBROUTINE writeba(FUNIT)
1008  IMPLICIT NONE
1009 
1010 ! File unit
1011  INTEGER, OPTIONAL, INTENT(IN) :: FUNIT
1012 
1013  INTEGER M, N
1014 
1015 ! Loop counter for continuum and discrete species
1016  INTEGER Nsp
1017 
1018  IF(.NOT.PRESENT(funit)) THEN
1019  write(*,1000)
1020  ELSE
1021  write(funit,1000)
1022  ENDIF
1023 
1024 ! Initialize counters
1025  nsp = 0
1026 
1027  DO n = 1, lng
1028  nsp = nsp + 1
1029  IF(lba(nsp)) THEN
1030  IF(.NOT.PRESENT(funit)) THEN
1031  write(*,1001)n
1032  ELSE
1033  write(funit,1001) n
1034  ENDIF
1035  ENDIF
1036  ENDDO
1037 
1038 ! Compare aliaes between solids phases
1039  DO m = 1, lm
1040  DO n = 1, lns(m)
1041  nsp = nsp + 1
1042 
1043  IF(lba(nsp)) THEN
1044  IF(.NOT.PRESENT(funit)) THEN
1045  write(*,1002)m, n
1046  ELSE
1047  write(funit,1002)m, n
1048  ENDIF
1049  ENDIF
1050  ENDDO
1051  ENDDO
1052 
1053  RETURN
1054 
1055  1000 FORMAT(' Species aliases were not provided for the following:')
1056  1001 FORMAT(3x, ' Gas phase species ',i2)
1057  1002 FORMAT(3x, ' Solid phase ',i2,' specie ',i2)
1058 
1059  END SUBROUTINE writeba
1060 
1061  END SUBROUTINE mapaliases
1062 
1063 
1064 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
1065 ! Function name: setFlag !
1066 ! !
1067 ! Purpose: Creates a flag pointing to a particular string location. !
1068 ! !
1069 ! Variables referenced: None !
1070 ! !
1071 !~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~!
1072  CHARACTER(len=512) FUNCTION setflag(fill, flg1, flg2) RESULT(OUT)
1074  IMPLICIT NONE
1075 
1076 ! Pass Arguments:
1077 !---------------------------------------------------------------------//
1078 ! Number of leading spaces to fill with dashes. This is used to jump
1079 ! past any lead-in text.
1080  INTEGER, INTENT(IN) :: fill
1081 ! Location in string to place the pointer.
1082  INTEGER, INTENT(IN) :: flg1
1083 ! Optional - location of second pointer in string
1084  INTEGER, INTENT(IN), OPTIONAL :: flg2
1085 
1086 ! Local Variables:
1087 !---------------------------------------------------------------------//
1088  INTEGER L, FILL1, FILL2
1089 
1090 ! Create a string with "FILL" dash characters.
1091  out = ''
1092  DO l = 1, fill-1
1093  WRITE(out,"(A,A)") trim(out), '-'
1094  ENDDO
1095 
1096 ! If a second pointer is present, determined the larger of the two.
1097  IF(PRESENT(flg2)) THEN
1098  IF(flg1 < flg2) THEN
1099  fill1 = flg1 - 1
1100  fill2 = (flg2-flg1) - 1
1101  ELSE
1102  fill1 = flg2 - 1
1103  fill2 = (flg1-flg2) - 1
1104  ENDIF
1105  ELSE
1106  fill1 = flg1 - 1
1107  fill2 = 0
1108  ENDIF
1109 
1110 ! Fill with dashes up to the the first pointer. ----^
1111  DO l = 1, fill1
1112  WRITE(out,"(A,A)") trim(out), '-'
1113  ENDDO
1114  WRITE(out,"(A,A)") trim(out), '^'
1115 ! Fill with dashes up to the second pointer. ----^---^
1116  IF(fill2 > 0) THEN
1117  DO l = 1, fill2
1118  WRITE(out,"(A,A)") trim(out), '-'
1119  ENDDO
1120  WRITE(out,"(A,A)") trim(out), '^'
1121  ENDIF
1122 
1123  END FUNCTION setflag
1124 
1125 END MODULE parse
subroutine writeba(FUNIT)
Definition: parse_mod.f:1007
integer, parameter dim_n_g
Definition: param_mod.f:69
subroutine setreaction(RxN, lNg, lSAg, lM, lNs, lSAs, lDH, lfDH)
Definition: parse_mod.f:71
subroutine splitaliasandcoeff(lName, lChemEq, lStart, lEnd, lAlias, lCoeff)
Definition: parse_mod.f:578
character(len=32), dimension(:), allocatable rxn_name
Definition: parse_mod.f:34
logical function compare(V1, V2)
Definition: toleranc_mod.f:94
double precision, parameter one
Definition: param1_mod.f:29
subroutine splitentries(lName, lChemEq, lStart, lEnd, lNo, lAlias, lCoeff)
Definition: parse_mod.f:494
logical des_rxn
Definition: parse_mod.f:23
logical more_chemeq
Definition: parse_mod.f:31
character(len=3), parameter end_blk
Definition: parse_mod.f:18
logical function checkmatch(lSA, ceSA)
Definition: parse_mod.f:936
double precision, dimension(:,:), allocatable des_usrfdh
Definition: parse_mod.f:54
integer, parameter dim_m
Definition: param_mod.f:67
character(len=8), parameter des_rxn_blk
Definition: parse_mod.f:17
subroutine updatemap(lnP, lpMap, llNoP)
Definition: parse_mod.f:970
double precision, parameter undefined
Definition: param1_mod.f:18
character(len=4), parameter rxn_blk
Definition: parse_mod.f:16
logical in_construct
Definition: parse_mod.f:28
subroutine checksplit(lName, lChemEq, lrEnd, lpStart, lSkip)
Definition: parse_mod.f:287
character(len=512) function setflag(fill, flg1, flg2)
Definition: parse_mod.f:1073
character(len=32), dimension(:), allocatable des_rxn_name
Definition: parse_mod.f:48
integer pe_io
Definition: compar_mod.f:30
subroutine mapaliases(lName, lChemEq, lNg, lSAg, lM, lNs, lSAs, lNo, lAlias, lCoeff, lSgn, lStart, lBA, lRxN)
Definition: parse_mod.f:816
subroutine mfix_exit(myID, normal_termination)
Definition: exit.f:5
subroutine checkblankaliases(lNg, lSAg, lM, lNs, lSAs, lBA)
Definition: parse_mod.f:750
Definition: exit.f:2
double precision, dimension(:), allocatable usrdh
Definition: parse_mod.f:38
double precision, dimension(:), allocatable des_usrdh
Definition: parse_mod.f:52
integer, parameter unit_log
Definition: funits_mod.f:21
Definition: param_mod.f:2
Definition: parse_mod.f:1
character(len=512), dimension(:), allocatable des_rxn_chem_eq
Definition: parse_mod.f:50
character(len=1), parameter end_str
Definition: parse_mod.f:13
double precision, dimension(:,:), allocatable usrfdh
Definition: parse_mod.f:40
integer mype
Definition: compar_mod.f:24
character(len=512), dimension(:), allocatable rxn_chem_eq
Definition: parse_mod.f:36
integer, parameter undefined_i
Definition: param1_mod.f:19
logical reading_rxn
Definition: parse_mod.f:20
integer, parameter dim_n_s
Definition: param_mod.f:71
logical reading_rate
Definition: parse_mod.f:21
logical in_des_construct
Definition: parse_mod.f:45
character(len=2), parameter start_str
Definition: parse_mod.f:12
logical tfm_rxn
Definition: parse_mod.f:24
subroutine make_upper_case(LINE_STRING, MAXCOL)
double precision, parameter zero
Definition: param1_mod.f:27