Line data Source code
1 : ! Include shortname defintions, so that the F77 code does not have to be modified to
2 : ! reference the CARMA structure.
3 : #include "carma_globaer.h"
4 :
5 : !! This routine checks if the coremass exceeds the total.
6 : !!
7 : !! NOTE: The fixer in this case is not mass conserving, and thus can have
8 : !! the effect of creating mass of the concentration element (by removing
9 : !! the negative values).
10 : !!
11 : !! NOTE: Errors will only be logged and runs aborted if they are greater
12 : !! than roundoff error.
13 : !!
14 : !! @author Yunqian Zhu, Charles Bardeen
15 : !! @version Apr-2021
16 4169515397 : subroutine coremasscheck(carma, cstate, iz, fixcoremass,logmsg,abort, packagename, rc)
17 :
18 : ! types
19 : use carma_precision_mod
20 : use carma_enums_mod
21 : use carma_constants_mod
22 : use carma_types_mod
23 : use carmastate_mod
24 : use carma_mod
25 :
26 : implicit none
27 :
28 : type(carma_type), intent(inout) :: carma !! the carma object
29 : type(carmastate_type), intent(inout) :: cstate !! the carma state object
30 : integer, intent(in) :: iz
31 : logical, intent(in) :: fixcoremass !! if we fix the coremass if exceeds total
32 : logical, intent(in) :: logmsg !! write message to log file if core mass check fails
33 : logical, intent(in) :: abort !! set return code to ERROR and return if if core mass check fails
34 : character(len=*), intent(in) :: packagename !! string to add to error message
35 : integer, intent(inout) :: rc !! return code, negative
36 :
37 : ! Local declarations
38 : integer :: igroup,ibin
39 : integer :: iepart,i,icore
40 : real(kind=f) :: coremass
41 : real(kind=f) :: factor
42 : real(kind=f),parameter :: roundoff = 1.e-14_f
43 :
44 : ! check the coremass exceeding the total mass
45 12508546191 : do igroup = 1,NGROUP
46 :
47 : ! Only check groups that have coremasses
48 12508546191 : if (ncore(igroup) > 0) then
49 :
50 4169515397 : iepart = ienconc(igroup)
51 :
52 87559823337 : do ibin = 1,NBIN
53 87559823337 : if (pc(iz, ibin, iepart) > 0._f) then
54 :
55 83390307940 : coremass = 0._f
56 :
57 >50034*10^7 : do i = 1, ncore(igroup)
58 >41695*10^7 : icore = icorelem(i,igroup)
59 >50034*10^7 : coremass = coremass + pc(iz, ibin, icore)
60 : end do ! i = 1, ncore(igroup)
61 :
62 83390307940 : if (coremass > pc(iz, ibin, iepart) * rmass(ibin, igroup)) then
63 83371819 : if (((coremass - pc(iz, ibin, iepart) * rmass(ibin, igroup)) / coremass) .gt. roundoff) then
64 5401109 : if (logmsg) then
65 0 : write(LUNOPRT,*) "Error - coremass exceeds total: ",packagename
66 0 : write(LUNOPRT,*) "coremass",coremass,"total",pc(iz, ibin, iepart) * rmass(ibin, igroup)
67 : end if
68 :
69 5401109 : if (abort) then
70 0 : rc = RC_ERROR
71 : return
72 : end if
73 :
74 : ! Only fix large errors if requested.
75 5401109 : if (fixcoremass) then
76 5401109 : pc(iz, ibin, iepart) = coremass / rmass(ibin, igroup)
77 : endif
78 :
79 : ! Automatically fix roundoff sized errors, regardless of settings.
80 : else
81 77970710 : pc(iz, ibin, iepart) = coremass / rmass(ibin, igroup)
82 : end if
83 : end if
84 : end if
85 : end do
86 : end if
87 : end do
88 :
89 4169515397 : end subroutine coremasscheck
|