1 | ! |
---|
2 | module module_fr_sfire_phys |
---|
3 | |
---|
4 | use module_model_constants, only: cp,xlv |
---|
5 | use module_fr_sfire_util |
---|
6 | |
---|
7 | PRIVATE |
---|
8 | |
---|
9 | ! subroutines and functions |
---|
10 | PUBLIC:: init_fuel_cats,fire_ros,heat_fluxes,set_nfuel_cat,set_fire_params,write_fuels_m |
---|
11 | |
---|
12 | PUBLIC::fire_params |
---|
13 | |
---|
14 | ! arrays passed to fire_ros |
---|
15 | type fire_params |
---|
16 | real,pointer,dimension(:,:):: vx,vy ! wind velocity (m/s) |
---|
17 | real,pointer,dimension(:,:):: zsf ! terrain height (m) |
---|
18 | real,pointer,dimension(:,:):: dzdxf,dzdyf ! terrain grad (1) |
---|
19 | real,pointer,dimension(:,:):: bbb,betafl,phiwc,r_0 ! spread formula coefficients |
---|
20 | real,pointer,dimension(:,:):: fgip ! init mass of surface fuel (kg/m^2) |
---|
21 | real,pointer,dimension(:,:):: ischap ! if fuel is chaparral and want rate of spread treated differently |
---|
22 | end type fire_params |
---|
23 | |
---|
24 | ! use as |
---|
25 | ! type(fire_params)::fp |
---|
26 | |
---|
27 | !D in col 2 means quantity derived from the others |
---|
28 | ! |
---|
29 | ! Scalar constants (data same for all fuel categories): |
---|
30 | ! HFGL SURFACE FIRE HEAT FLUX THRESHOLD TO IGNITE CANOPY (W/m^2) |
---|
31 | ! CMBCNST JOULES PER KG OF DRY FUEL |
---|
32 | ! FUELHEAT FUEL PARTICLE LOW HEAT CONTENT, BTU/LB |
---|
33 | ! FUELMC_G FUEL PARTICLE (SURFACE) MOISTURE CONTENT |
---|
34 | !D BMST RATIO OF LATENT TO SENSIBLE HEAT FROM SFC BURN: |
---|
35 | ! % of total fuel mass that is water (not quite |
---|
36 | ! = % fuel moisture). BMST= (H20)/(H20+DRY) |
---|
37 | ! so BMST = FUELMC_G / (1 + FUELMC_G) where |
---|
38 | ! FUELMC_G = moisture content of surface fuel |
---|
39 | ! |
---|
40 | ! Data arrays indexed by fuel category: |
---|
41 | ! FGI INITIAL TOTAL MASS OF SURFACE FUEL (KG/M**2) |
---|
42 | ! FUELDEPTHM FUEL DEPTH, IN M (CONVERTED TO FT) |
---|
43 | ! SAVR FUEL PARTICLE SURFACE-AREA-TO-VOLUME RATIO, 1/FT |
---|
44 | ! GRASS: 3500., 10 hr fuel: 109., 100 hr fuel: 30. |
---|
45 | ! FUELMCE MOISTURE CONTENT OF EXTINCTION; 0.30 FOR MANY DEAD FUELS; 0.15 FOR GRASS |
---|
46 | ! FUELDENS OVENDRY PARTICLE DENSITY, LB/FT^3 |
---|
47 | ! ST FUEL PARTICLE TOTAL MINERAL CONTENT |
---|
48 | ! SE FUEL PARTICLE EFFECTIVE MINERAL CONTENT |
---|
49 | ! WEIGHT WEIGHTING PARAMETER THAT DETERMINES THE SLOPE OF THE MASS LOSS CURVE |
---|
50 | ! RANGES FROM ~5 (FAST BURNUP) TO 1000 ( ~40% DECR OVER 10 MIN). |
---|
51 | ! FCI_D INITIAL DRY MASS OF CANOPY FUEL |
---|
52 | ! FCT BURN OUT TIME FOR CANOPY FUEL, AFTER DRY (S) |
---|
53 | ! ichap Set=1 if fuel is chaparral and want the rate of spread treated differently, 0 if not |
---|
54 | !D FCI INITIAL TOTAL MASS OF CANOPY FUEL |
---|
55 | !D FCBR FUEL CANOPY BURN RATE (KG/M**2/S) |
---|
56 | |
---|
57 | ! ============================================================================= |
---|
58 | ! Anderson 13 surface fire fuel models, along with some |
---|
59 | ! estimated canopy properties (for crown fire). |
---|
60 | ! ============================================================================= |
---|
61 | ! --- Grass-dominated fuel models |
---|
62 | ! FUEL MODEL 1: Short grass (1 ft) |
---|
63 | ! FUEL MODEL 2: Timber (grass and understory) |
---|
64 | ! FUEL MODEL 3: Tall grass (2.5 ft) |
---|
65 | ! --- Shrub-dominated fuel models |
---|
66 | ! FUEL MODEL 4: Chaparral (6 ft) |
---|
67 | ! FUEL MODEL 5: Brush (2 ft) |
---|
68 | ! FUEL MODEL 6: Dormant brush, hardwood slash |
---|
69 | ! FUEL MODEL 7: Southern rough |
---|
70 | ! --- Forest litter-dominated fuel models |
---|
71 | ! FUEL MODEL 8: Closed timber litter |
---|
72 | ! FUEL MODEL 9: Hardwood litter |
---|
73 | ! FUEL MODEL 10: Timber (litter + understory) |
---|
74 | ! --- Logging debris-dominated fuel models |
---|
75 | ! FUEL MODEL 11: Light logging slash |
---|
76 | ! FUEL MODEL 12: Medium logging slash |
---|
77 | ! FUEL MODEL 13: Heavy logging slash |
---|
78 | ! --- Fuel-free |
---|
79 | ! FUEL MODEL 14: no fuel |
---|
80 | |
---|
81 | ! scalar fuel coefficients |
---|
82 | REAL, SAVE:: cmbcnst,hfgl,fuelmc_g,fuelmc_c |
---|
83 | ! computed values |
---|
84 | REAL, SAVE:: bmst,fuelheat |
---|
85 | |
---|
86 | ! defaults, may be changed in init_fuel_cats |
---|
87 | DATA cmbcnst / 17.433e+06/ ! J/kg dry fuel |
---|
88 | DATA hfgl / 17.e4 / ! W/m^2 |
---|
89 | DATA fuelmc_g / 0.08 / ! set = 0 for dry surface fuel |
---|
90 | DATA fuelmc_c / 1.00 / ! set = 0 for dry canopy |
---|
91 | ! REAL, PARAMETER :: bmst = fuelmc_g/(1+fuelmc_g) |
---|
92 | ! REAL, PARAMETER :: fuelheat = cmbcnst * 4.30e-04 ! convert J/kg to BTU/lb |
---|
93 | ! real, parameter :: xlv = 2.5e6 ! to make it selfcontained |
---|
94 | ! real, parameter :: cp = 7.*287./2 ! to make it selfcontained |
---|
95 | |
---|
96 | |
---|
97 | ! fuel categorytables |
---|
98 | INTEGER, PARAMETER :: nf=14 ! number of fuel categories in data stmts |
---|
99 | INTEGER, SAVE :: nfuelcats = 13 ! number of fuel categories that are specified |
---|
100 | INTEGER, PARAMETER :: mfuelcats = 30 ! allowable number of fuel categories |
---|
101 | INTEGER, PARAMETER :: zf = mfuelcats-nf ! number of zero fillers in data stmt |
---|
102 | INTEGER, SAVE :: no_fuel_cat = 14 ! special category outside of 1:nfuelcats |
---|
103 | CHARACTER (len=80), DIMENSION(mfuelcats ), save :: fuel_name |
---|
104 | INTEGER, DIMENSION( mfuelcats ), save :: ichap |
---|
105 | REAL , DIMENSION( mfuelcats ), save :: windrf,weight,fgi,fci,fci_d,fct,fcbr, & |
---|
106 | fueldepthm,fueldens,fuelmce, & |
---|
107 | savr,st,se |
---|
108 | DATA windrf /0.36, 0.36, 0.44, 0.55, 0.42, 0.44, 0.44, & |
---|
109 | 0.36, 0.36, 0.36, 0.36, 0.43, 0.46, 1e-7, zf*0 / |
---|
110 | DATA fgi / 0.166, 0.896, 0.674, 3.591, 0.784, 1.344, 1.091, & |
---|
111 | 1.120, 0.780, 2.692, 2.582, 7.749, 13.024, 1.e-7, zf*0. / |
---|
112 | DATA fueldepthm /0.305, 0.305, 0.762, 1.829, 0.61, 0.762,0.762, & |
---|
113 | 0.0610, 0.0610, 0.305, 0.305, 0.701, 0.914, 0.305,zf*0. / |
---|
114 | DATA savr / 3500., 2784., 1500., 1739., 1683., 1564., 1562., & |
---|
115 | 1889., 2484., 1764., 1182., 1145., 1159., 3500., zf*0. / |
---|
116 | DATA fuelmce / 0.12, 0.15, 0.25, 0.20, 0.20, 0.25, 0.40, & |
---|
117 | 0.30, 0.25, 0.25, 0.15, 0.20, 0.25, 0.12 , zf*0. / |
---|
118 | DATA fueldens / nf * 32., zf*0. / ! 32 if solid, 19 if rotten. |
---|
119 | DATA st / nf* 0.0555 , zf*0./ |
---|
120 | DATA se / nf* 0.010 , zf*0./ |
---|
121 | ! ----- Notes on weight: (4) - best fit of data from D. Latham (pers. comm.); |
---|
122 | ! (5)-(7) could be 60-120; (8)-(10) could be 300-1600; |
---|
123 | ! (11)-(13) could be 300-1600 |
---|
124 | DATA weight / 7., 7., 7., 180., 100., 100., 100., & |
---|
125 | 900., 900., 900., 900., 900., 900., 7. , zf*0./ |
---|
126 | ! ----- 1.12083 is 5 tons/acre. 5-50 tons/acre orig., 100-300 after blowdown |
---|
127 | DATA fci_d / 0., 0., 0., 1.123, 0., 0., 0., & |
---|
128 | 1.121, 1.121, 1.121, 1.121, 1.121, 1.121, 0., zf*0./ |
---|
129 | DATA fct / 60., 60., 60., 60., 60., 60., 60., & |
---|
130 | 60., 120., 180., 180., 180., 180. , 60. , zf*0. / |
---|
131 | DATA ichap / 0, 0, 0, 1, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0 , zf*0/ |
---|
132 | ! ========================================================================= |
---|
133 | |
---|
134 | contains |
---|
135 | |
---|
136 | subroutine init_fuel_cats |
---|
137 | implicit none |
---|
138 | !*** purpose: initialize fuel tables and variables by constants |
---|
139 | !*** arguments: none |
---|
140 | logical, external:: wrf_dm_on_monitor |
---|
141 | !$ integer, external:: OMP_GET_THREAD_NUM |
---|
142 | !*** local |
---|
143 | integer:: i,j,k,ii,iounit |
---|
144 | character(len=128):: msg |
---|
145 | !*** executable |
---|
146 | |
---|
147 | ! read |
---|
148 | namelist /fuel_scalars/ cmbcnst,hfgl,fuelmc_g,fuelmc_c,nfuelcats,no_fuel_cat |
---|
149 | namelist /fuel_categories/ fuel_name,windrf,fgi,fueldepthm,savr, & |
---|
150 | fuelmce,fueldens,st,se,weight,fci_d,fct,ichap |
---|
151 | |
---|
152 | !$ if (OMP_GET_THREAD_NUM() .ne. 0)then |
---|
153 | !$ call crash('init_fuel_cats: must be called from master thread') |
---|
154 | !$ endif |
---|
155 | |
---|
156 | IF ( wrf_dm_on_monitor() ) THEN |
---|
157 | ! if we are the master task, read the file |
---|
158 | iounit=open_text_file('namelist.fire','read') |
---|
159 | read(iounit,fuel_scalars) |
---|
160 | read(iounit,fuel_categories) |
---|
161 | CLOSE(iounit) |
---|
162 | |
---|
163 | if (nfuelcats>mfuelcats) then |
---|
164 | write(msg,*)'nfuelcats=',nfuelcats,' is too large, increase mfuelcats' |
---|
165 | call crash(msg) |
---|
166 | endif |
---|
167 | if (no_fuel_cat >= 1 .and. no_fuel_cat <= nfuelcats)then |
---|
168 | write(msg,*)'no_fuel_cat=',no_fuel_cat,' may not be between 1 and nfuelcats=',nfuelcats |
---|
169 | call crash(msg) |
---|
170 | endif |
---|
171 | ENDIF |
---|
172 | |
---|
173 | ! broadcast the contents of the file |
---|
174 | call wrf_dm_bcast_real(cmbcnst,1) |
---|
175 | call wrf_dm_bcast_real(hfgl,1) |
---|
176 | call wrf_dm_bcast_real(fuelmc_g,1) |
---|
177 | call wrf_dm_bcast_real(fuelmc_c,1) |
---|
178 | call wrf_dm_bcast_integer(nfuelcats,1) |
---|
179 | call wrf_dm_bcast_integer(no_fuel_cat,1) |
---|
180 | call wrf_dm_bcast_real(windrf, nfuelcats) |
---|
181 | call wrf_dm_bcast_real(fgi, nfuelcats) |
---|
182 | call wrf_dm_bcast_real(fueldepthm,nfuelcats) |
---|
183 | call wrf_dm_bcast_real(savr, nfuelcats) |
---|
184 | call wrf_dm_bcast_real(fuelmce, nfuelcats) |
---|
185 | call wrf_dm_bcast_real(fueldens, nfuelcats) |
---|
186 | call wrf_dm_bcast_real(st, nfuelcats) |
---|
187 | call wrf_dm_bcast_real(se, nfuelcats) |
---|
188 | call wrf_dm_bcast_real(weight, nfuelcats) |
---|
189 | call wrf_dm_bcast_real(fci_d, nfuelcats) |
---|
190 | call wrf_dm_bcast_real(fct, nfuelcats) |
---|
191 | call wrf_dm_bcast_integer(ichap, nfuelcats) |
---|
192 | |
---|
193 | ! compute derived scalars |
---|
194 | |
---|
195 | bmst = fuelmc_g/(1+fuelmc_g) |
---|
196 | fuelheat = cmbcnst * 4.30e-04 ! convert J/kg to BTU/lb |
---|
197 | |
---|
198 | ! compute derived fuel category coefficients |
---|
199 | |
---|
200 | DO i = 1,nfuelcats |
---|
201 | fci(i) = (1.+fuelmc_c)*fci_d(i) |
---|
202 | if(fct(i) .ne. 0.)then |
---|
203 | fcbr(i) = fci_d(i)/fct(i) ! avoid division by zero |
---|
204 | else |
---|
205 | fcbr(i) = 0 |
---|
206 | endif |
---|
207 | END DO |
---|
208 | |
---|
209 | ! prints |
---|
210 | |
---|
211 | call message('**********************************************************') |
---|
212 | call message('FUEL COEFFICIENTS') |
---|
213 | write(msg,8)'cmbcnst ',cmbcnst |
---|
214 | call message(msg) |
---|
215 | write(msg,8)'hfgl ',hfgl |
---|
216 | call message(msg) |
---|
217 | write(msg,8)'fuelmc_g ',fuelmc_g |
---|
218 | call message(msg) |
---|
219 | write(msg,8)'fuelmc_c ',fuelmc_c |
---|
220 | call message(msg) |
---|
221 | write(msg,8)'bmst ',bmst |
---|
222 | call message(msg) |
---|
223 | write(msg,8)'fuelheat ',fuelheat |
---|
224 | call message(msg) |
---|
225 | write(msg,7)'nfuelcats ',nfuelcats |
---|
226 | call message(msg) |
---|
227 | write(msg,7)'no_fuel_cat',no_fuel_cat |
---|
228 | call message(msg) |
---|
229 | |
---|
230 | j=1 |
---|
231 | 7 format(a,5(1x,i8,4x)) |
---|
232 | 8 format(a,5(1x,g12.5e2)) |
---|
233 | 9 format(a,5(1x,a)) |
---|
234 | do i=1,nfuelcats,j |
---|
235 | k=min(i+j-1,nfuelcats) |
---|
236 | call message(' ') |
---|
237 | write(msg,7)'CATEGORY ',(ii,ii=i,k) |
---|
238 | call message(msg) |
---|
239 | write(msg,9)'fuel name ',(fuel_name(ii),ii=i,k) |
---|
240 | call message(msg) |
---|
241 | ! write(msg,8)'windrf ',(windrf(ii),ii=i,k) |
---|
242 | ! call message(msg) |
---|
243 | write(msg,8)'fgi ',(fgi(ii),ii=i,k) |
---|
244 | call message(msg) |
---|
245 | write(msg,8)'fueldepthm',(fueldepthm(ii),ii=i,k) |
---|
246 | call message(msg) |
---|
247 | write(msg,8)'savr ',(savr(ii),ii=i,k) |
---|
248 | call message(msg) |
---|
249 | write(msg,8)'fuelmce ',(fuelmce(ii),ii=i,k) |
---|
250 | call message(msg) |
---|
251 | write(msg,8)'fueldens ',(fueldens(ii),ii=i,k) |
---|
252 | call message(msg) |
---|
253 | write(msg,8)'st ',(st(ii),ii=i,k) |
---|
254 | call message(msg) |
---|
255 | write(msg,8)'se ',(se(ii),ii=i,k) |
---|
256 | call message(msg) |
---|
257 | write(msg,8)'weight ',(weight(ii),ii=i,k) |
---|
258 | call message(msg) |
---|
259 | write(msg,8)'fci_d ',(fci_d(ii),ii=i,k) |
---|
260 | call message(msg) |
---|
261 | write(msg,8)'fct ',(fct(ii),ii=i,k) |
---|
262 | call message(msg) |
---|
263 | write(msg,7)'ichap ',(ichap(ii),ii=i,k) |
---|
264 | call message(msg) |
---|
265 | write(msg,8)'fci ',(fci(ii),ii=i,k) |
---|
266 | call message(msg) |
---|
267 | write(msg,8)'fcbr ',(fcbr(ii),ii=i,k) |
---|
268 | call message(msg) |
---|
269 | enddo |
---|
270 | call message('**********************************************************') |
---|
271 | |
---|
272 | ! and print to file |
---|
273 | IF ( wrf_dm_on_monitor() ) THEN |
---|
274 | call write_fuels_m(61,30.,1.) |
---|
275 | ENDIF |
---|
276 | end subroutine init_fuel_cats |
---|
277 | |
---|
278 | |
---|
279 | subroutine write_fuels_m(nsteps,maxwind,maxslope) |
---|
280 | implicit none |
---|
281 | integer, intent(in):: nsteps ! number of steps for speed computation |
---|
282 | real, intent(in):: maxwind,maxslope ! computer from zero to these |
---|
283 | |
---|
284 | integer:: iounit,k,j,i |
---|
285 | type(fire_params)::fp |
---|
286 | !type fire_params |
---|
287 | !real,pointer,dimension(:,:):: vx,vy ! wind velocity (m/s) |
---|
288 | !real,pointer,dimension(:,:):: zsf ! terrain height (m) |
---|
289 | !real,pointer,dimension(:,:):: dzdxf,dzdyf ! terrain grad (1) |
---|
290 | !real,pointer,dimension(:,:):: bbb,betafl,phiwc,r_0 ! spread formula coefficients |
---|
291 | !real,pointer,dimension(:,:):: fgip ! init mass of surface fuel (kg/m^2) |
---|
292 | !real,pointer,dimension(:,:):: ischap ! 1 if chapparal |
---|
293 | !end type fire_params |
---|
294 | real, dimension(1:2,1:nsteps), target::vx,vy,zsf,dzdxf,dzdyf,bbb,betafl,phiwc,r_0,fgip,ischap |
---|
295 | real, dimension(1:2,1:nsteps)::nfuel_cat,fuel_time,ros |
---|
296 | real::ros_base,ros_wind,ros_slope,propx,propy,r |
---|
297 | |
---|
298 | fp%vx=>vx |
---|
299 | fp%vy=>vy |
---|
300 | fp%dzdxf=>dzdxf |
---|
301 | fp%dzdyf=>dzdyf |
---|
302 | fp%bbb=>bbb |
---|
303 | fp%betafl=>betafl |
---|
304 | fp%phiwc=>phiwc |
---|
305 | fp%r_0=>r_0 |
---|
306 | fp%fgip=>fgip |
---|
307 | fp%ischap=>ischap |
---|
308 | |
---|
309 | iounit = open_text_file('fuels.m','write') |
---|
310 | |
---|
311 | 10 format('fuel(',i3,').',a,'=',"'",a,"'",';% ',a) |
---|
312 | do k=1,nfuelcats |
---|
313 | write(iounit,10)k,'fuel_name',trim(fuel_name(k)),'FUEL MODEL NAME' |
---|
314 | call write_var(k,'windrf',windrf(k),'WIND REDUCTION FACTOR FROM 20ft TO MIDFLAME HEIGHT' ) |
---|
315 | call write_var(k,'fgi',fgi(k),'INITIAL TOTAL MASS OF SURFACE FUEL (KG/M**2)' ) |
---|
316 | call write_var(k,'fueldepthm',fueldepthm(k),'FUEL DEPTH (M)') |
---|
317 | call write_var(k,'savr',savr(k),'FUEL PARTICLE SURFACE-AREA-TO-VOLUME RATIO, 1/FT') |
---|
318 | call write_var(k,'fuelmce',fuelmce(k),'MOISTURE CONTENT OF EXTINCTION') |
---|
319 | call write_var(k,'fueldens',fueldens(k),'OVENDRY PARTICLE DENSITY, LB/FT^3') |
---|
320 | call write_var(k,'st',st(k),'FUEL PARTICLE TOTAL MINERAL CONTENT') |
---|
321 | call write_var(k,'se',se(k),'FUEL PARTICLE EFFECTIVE MINERAL CONTENT') |
---|
322 | call write_var(k,'weight',weight(k),'WEIGHTING PARAMETER THAT DETERMINES THE SLOPE OF THE MASS LOSS CURVE') |
---|
323 | call write_var(k,'fci_d',fci_d(k),'INITIAL DRY MASS OF CANOPY FUEL') |
---|
324 | call write_var(k,'fct',fct(k),'BURN OUT TIME FOR CANOPY FUEL, AFTER DRY (S)') |
---|
325 | call write_var(k,'ichap',float(ichap(k)),'1 if chaparral, 0 if not') |
---|
326 | call write_var(k,'fci',fci(k),'INITIAL TOTAL MASS OF CANOPY FUEL') |
---|
327 | call write_var(k,'fcbr',fcbr(k),'FUEL CANOPY BURN RATE (KG/M**2/S)') |
---|
328 | call write_var(k,'hfgl',hfgl,'SURFACE FIRE HEAT FLUX THRESHOLD TO IGNITE CANOPY (W/m^2)') |
---|
329 | call write_var(k,'cmbcnst',cmbcnst,'JOULES PER KG OF DRY FUEL') |
---|
330 | call write_var(k,'fuelheat',fuelheat,'FUEL PARTICLE LOW HEAT CONTENT, BTU/LB') |
---|
331 | call write_var(k,'fuelmc_g',fuelmc_g,'FUEL PARTICLE (SURFACE) MOISTURE CONTENT') |
---|
332 | call write_var(k,'fuelmc_c',fuelmc_c,'FUEL PARTICLE (CANOPY) MOISTURE CONTENT') |
---|
333 | ! set up fuel arrays |
---|
334 | !subroutine set_fire_params( & |
---|
335 | ! ifds,ifde,jfds,jfde, & |
---|
336 | ! ifms,ifme,jfms,jfme, & |
---|
337 | ! ifts,ifte,jfts,jfte, & |
---|
338 | ! fdx,fdy,nfuel_cat0, & |
---|
339 | ! nfuel_cat,fuel_time, & |
---|
340 | ! fp ) |
---|
341 | nfuel_cat = k |
---|
342 | call set_fire_params( & |
---|
343 | 1,2,1,nsteps, & |
---|
344 | 1,2,1,nsteps, & |
---|
345 | 1,2,1,nsteps, & |
---|
346 | 0.,0.,k, & |
---|
347 | nfuel_cat,fuel_time, & |
---|
348 | fp ) |
---|
349 | ! set up windspeed and slope table |
---|
350 | propx=1. |
---|
351 | propy=0. |
---|
352 | do j=1,nsteps |
---|
353 | r=float(j-1)/float(nsteps-1) |
---|
354 | ! line 1 varies windspeed (in x direction), zero slope |
---|
355 | vx(1,j)=maxwind*r |
---|
356 | vy(1,j)=0. |
---|
357 | dzdxf(1,j)=0. |
---|
358 | dzdyf(1,j)=0. |
---|
359 | ! line 2 varies slope (in x direction), zero slope |
---|
360 | vx(2,j)=0. |
---|
361 | vy(2,j)=0. |
---|
362 | dzdxf(2,j)=maxslope*r |
---|
363 | dzdyf(2,j)=0. |
---|
364 | enddo |
---|
365 | do j=1,nsteps |
---|
366 | do i=1,2 |
---|
367 | call fire_ros(ros_base,ros_wind,ros_slope, & |
---|
368 | propx,propy,i,j,fp) |
---|
369 | ros(i,j)=ros_base+ros_wind+ros_slope |
---|
370 | enddo |
---|
371 | write(iounit,13)k,'wind',j,vx(1,j),'wind speed' |
---|
372 | write(iounit,13)k,'ros_wind',j,ros(1,j),'rate of spread for the wind speed' |
---|
373 | write(iounit,13)k,'slope',j,dzdxf(2,j),'slope' |
---|
374 | write(iounit,13)k,'ros_slope',j,ros(2,j),'rate of spread for the slope' |
---|
375 | enddo |
---|
376 | enddo |
---|
377 | 13 format('fuel(',i3,').',a,'(',i3,')=',g12.5e2,';% ',a) |
---|
378 | |
---|
379 | close(iounit) |
---|
380 | ! stop |
---|
381 | |
---|
382 | contains |
---|
383 | |
---|
384 | subroutine write_var(k,name,value,descr) |
---|
385 | ! write entry for one variable |
---|
386 | integer, intent(in)::k |
---|
387 | character(len=*), intent(in)::name,descr |
---|
388 | real, intent(in)::value |
---|
389 | write(iounit,11)k,name,value |
---|
390 | write(iounit,12)k,name,descr |
---|
391 | 11 format('fuel(',i3,').',a,'=',g12.5e2, ';') |
---|
392 | 12 format('fuel(',i3,').',a,"_descr='",a,"';") |
---|
393 | end subroutine write_var |
---|
394 | |
---|
395 | end subroutine write_fuels_m |
---|
396 | |
---|
397 | ! |
---|
398 | !******************* |
---|
399 | ! |
---|
400 | |
---|
401 | subroutine set_fire_params( & |
---|
402 | ifds,ifde,jfds,jfde, & |
---|
403 | ifms,ifme,jfms,jfme, & |
---|
404 | ifts,ifte,jfts,jfte, & |
---|
405 | fdx,fdy,nfuel_cat0, & |
---|
406 | nfuel_cat,fuel_time, & |
---|
407 | fp ) |
---|
408 | |
---|
409 | implicit none |
---|
410 | |
---|
411 | !*** purpose: Set all fire model params arrays, constant values. |
---|
412 | |
---|
413 | !*** arguments |
---|
414 | integer, intent(in)::ifds,ifde,jfds,jfde ! fire domain bounds |
---|
415 | integer, intent(in)::ifts,ifte,jfts,jfte ! fire tile bounds |
---|
416 | integer, intent(in)::ifms,ifme,jfms,jfme ! memory array bounds |
---|
417 | real, intent(in):: fdx,fdy ! fire mesh spacing |
---|
418 | integer,intent(in)::nfuel_cat0 ! default fuel category, if nfuel_cat=0 |
---|
419 | real, intent(in),dimension(ifms:ifme, jfms:jfme)::nfuel_cat ! fuel data |
---|
420 | real, intent(out), dimension(ifms:ifme, jfms:jfme)::fuel_time ! fire params arrays |
---|
421 | type(fire_params),intent(inout)::fp |
---|
422 | |
---|
423 | !*** local |
---|
424 | |
---|
425 | real:: fuelload, fueldepth, rtemp1, rtemp2, & |
---|
426 | qig, epsilon, rhob, wn, betaop, e, c, & |
---|
427 | xifr, etas, etam, a, gammax, gamma, ratio, ir, & |
---|
428 | fuelloadm,fdxinv,fdyinv |
---|
429 | integer:: i,j,k |
---|
430 | integer::nerr |
---|
431 | character(len=128)::msg |
---|
432 | |
---|
433 | !*** executable |
---|
434 | |
---|
435 | nerr=0 |
---|
436 | do j=jfts,jfte |
---|
437 | do i=ifts,ifte |
---|
438 | ! fuel category |
---|
439 | k=int( nfuel_cat(i,j) ) |
---|
440 | if(k.eq.no_fuel_cat)then ! no fuel |
---|
441 | fp%fgip(i,j)=0. ! no mass |
---|
442 | fp%ischap(i,j)=0. |
---|
443 | fp%betafl(i,j)=0. ! to prevent division by zero |
---|
444 | fp%bbb(i,j)=1. ! |
---|
445 | fuel_time(i,j)=7./0.85 ! does not matter, just what was there before |
---|
446 | fp%phiwc(i,j)=0. |
---|
447 | fp%r_0(i,j)=0. ! no fuel, no spread. |
---|
448 | else |
---|
449 | if(k.eq.0.and.nfuel_cat0.ge.1.and.nfuel_cat0.le.nfuelcats)then |
---|
450 | ! replace k=0 by default |
---|
451 | k=nfuel_cat0 |
---|
452 | nerr=nerr+1 |
---|
453 | endif |
---|
454 | |
---|
455 | if(k.lt.1.or.k.gt.nfuelcats)then |
---|
456 | !$OMP CRITICAL(SFIRE_PHYS_CRIT) |
---|
457 | write(msg,'(3(a,i5))')'nfuel_cat(', i ,',',j,')=',k |
---|
458 | !$OMP END CRITICAL(SFIRE_PHYS_CRIT) |
---|
459 | call message(msg) |
---|
460 | call crash('set_fire_params: fuel category out of bounds') |
---|
461 | endif |
---|
462 | |
---|
463 | fuel_time(i,j)=weight(k)/0.85 ! cell based |
---|
464 | |
---|
465 | ! set fuel time constant: weight=1000 => 40% decrease over 10 min |
---|
466 | ! fuel decreases as exp(-t/fuel_time) |
---|
467 | ! exp(-600*0.85/1000) = approx 0.6 |
---|
468 | |
---|
469 | fp%ischap(i,j)=ichap(k) |
---|
470 | fp%fgip(i,j)=fgi(k) |
---|
471 | |
---|
472 | ! ...Settings of fire spread parameters from Rothermel follows. These |
---|
473 | ! don't need to be recalculated later. |
---|
474 | |
---|
475 | fuelloadm= (1.-bmst) * fgi(k) ! fuelload without moisture |
---|
476 | fuelload = fuelloadm * (.3048)**2 * 2.205 ! to lb/ft^2 |
---|
477 | fueldepth = fueldepthm(k)/0.3048 ! to ft |
---|
478 | fp%betafl(i,j) = fuelload/(fueldepth * fueldens(k))! packing ratio |
---|
479 | betaop = 3.348 * savr(k)**(-0.8189) ! optimum packing ratio |
---|
480 | qig = 250. + 1116.*fuelmc_g ! heat of preignition, btu/lb |
---|
481 | epsilon = exp(-138./savr(k) ) ! effective heating number |
---|
482 | rhob = fuelload/fueldepth ! ovendry bulk density, lb/ft^3 |
---|
483 | |
---|
484 | c = 7.47 * exp( -0.133 * savr(k)**0.55) ! const in wind coef |
---|
485 | fp%bbb(i,j) = 0.02526 * savr(k)**0.54 ! const in wind coef |
---|
486 | e = 0.715 * exp( -3.59e-4 * savr(k)) ! const in wind coef |
---|
487 | fp%phiwc(i,j) = c * (fp%betafl(i,j)/betaop)**(-e) |
---|
488 | |
---|
489 | rtemp2 = savr(k)**1.5 |
---|
490 | gammax = rtemp2/(495. + 0.0594*rtemp2) ! maximum rxn vel, 1/min |
---|
491 | a = 1./(4.774 * savr(k)**0.1 - 7.27) ! coef for optimum rxn vel |
---|
492 | ratio = fp%betafl(i,j)/betaop |
---|
493 | gamma = gammax *(ratio**a) *exp(a*(1.-ratio)) !optimum rxn vel, 1/min |
---|
494 | |
---|
495 | wn = fuelload/(1 + st(k)) ! net fuel loading, lb/ft^2 |
---|
496 | rtemp1 = fuelmc_g/fuelmce(k) |
---|
497 | etam = 1.-2.59*rtemp1 +5.11*rtemp1**2 -3.52*rtemp1**3 !moist damp coef |
---|
498 | etas = 0.174* se(k)**(-0.19) ! mineral damping coef |
---|
499 | ir = gamma * wn * fuelheat * etam * etas !rxn intensity,btu/ft^2 min |
---|
500 | ! irm = ir * 1055./( 0.3048**2 * 60.) * 1.e-6 !for mw/m^2 |
---|
501 | |
---|
502 | xifr = exp( (0.792 + 0.681*savr(k)**0.5) & |
---|
503 | * (fp%betafl(i,j)+0.1)) /(192. + 0.2595*savr(k)) ! propagating flux ratio |
---|
504 | |
---|
505 | ! ... r_0 is the spread rate for a fire on flat ground with no wind. |
---|
506 | fp%r_0(i,j) = ir*xifr/(rhob * epsilon *qig) ! default spread rate in ft/min |
---|
507 | |
---|
508 | endif |
---|
509 | enddo |
---|
510 | enddo |
---|
511 | |
---|
512 | if(nerr.gt.1)then |
---|
513 | !$OMP CRITICAL(SFIRE_PHYS_CRIT) |
---|
514 | write(msg,'(a,i6,a)')'set_fire_params: WARNING: fuel category 0 replaced in',nerr,' cells' |
---|
515 | !$OMP END CRITICAL(SFIRE_PHYS_CRIT) |
---|
516 | call message(msg) |
---|
517 | endif |
---|
518 | |
---|
519 | end subroutine set_fire_params |
---|
520 | |
---|
521 | ! |
---|
522 | !******************* |
---|
523 | ! |
---|
524 | |
---|
525 | subroutine heat_fluxes(dt, & |
---|
526 | ifms,ifme,jfms,jfme, & ! memory dims |
---|
527 | ifts,ifte,jfts,jfte, & ! tile dims |
---|
528 | iffs,iffe,jffs,jffe, & ! fuel_frac_burnt dims |
---|
529 | fgip,fuel_frac_burnt, & !in |
---|
530 | grnhft,grnqft) !out |
---|
531 | implicit none |
---|
532 | |
---|
533 | !*** purpose |
---|
534 | ! compute the heat fluxes on the fire grid cells |
---|
535 | |
---|
536 | !*** arguments |
---|
537 | real, intent(in)::dt ! dt the fire time step (the fire model advances time by this) |
---|
538 | integer, intent(in)::ifts,ifte,jfts,jfte,ifms,ifme,jfms,jfme,iffs,iffe,jffs,jffe ! dimensions |
---|
539 | real, intent(in),dimension(ifms:ifme,jfms:jfme):: fgip |
---|
540 | real, intent(in),dimension(iffs:iffe,jffs:jffe):: fuel_frac_burnt |
---|
541 | real, intent(out),dimension(ifms:ifme,jfms:jfme):: grnhft,grnqft |
---|
542 | |
---|
543 | !*** local |
---|
544 | integer::i,j |
---|
545 | real:: dmass |
---|
546 | |
---|
547 | !*** executable |
---|
548 | do j=jfts,jfte |
---|
549 | do i=ifts,ifte |
---|
550 | dmass = & ! surface fuel dry mass burnt this call (kg/m^2) |
---|
551 | fgip(i,j) & ! init mass from fuel model no (kg/m^2) = fgi(nfuel_cat(i,j) |
---|
552 | * fuel_frac_burnt(i,j) ! fraction burned this call (1) |
---|
553 | grnhft(i,j) = (dmass/dt)*(1.-bmst)*cmbcnst ! surface fire sensible heat flux W/m^2 |
---|
554 | grnqft(i,j) = (bmst+(1.-bmst)*.56)*(dmass/dt)*xlv ! surface fire latent heat flux W/m |
---|
555 | ! xlv is defined in module_model_constants.. Assume 56% of cellulose molecule mass is water. |
---|
556 | enddo |
---|
557 | enddo |
---|
558 | |
---|
559 | end subroutine heat_fluxes |
---|
560 | |
---|
561 | ! |
---|
562 | !********************** |
---|
563 | ! |
---|
564 | |
---|
565 | |
---|
566 | subroutine set_nfuel_cat( & |
---|
567 | ifms,ifme,jfms,jfme, & |
---|
568 | ifts,ifte,jfts,jfte, & |
---|
569 | ifuelread,nfuel_cat0,zsf,nfuel_cat) |
---|
570 | |
---|
571 | implicit none |
---|
572 | |
---|
573 | ! set fuel distributions for testing |
---|
574 | integer, intent(in):: ifts,ifte,jfts,jfte, & |
---|
575 | ifms,ifme,jfms,jfme |
---|
576 | |
---|
577 | integer, intent(in)::ifuelread,nfuel_cat0 |
---|
578 | real, intent(in), dimension(ifms:ifme, jfms:jfme)::zsf |
---|
579 | real, intent(out), dimension(ifms:ifme, jfms:jfme)::nfuel_cat |
---|
580 | |
---|
581 | !*** local |
---|
582 | |
---|
583 | ! parameters to control execution |
---|
584 | integer:: i,j,iu1 |
---|
585 | real:: t1 |
---|
586 | character(len=128)msg |
---|
587 | |
---|
588 | !$OMP CRITICAL(SFIRE_PHYS_CRIT) |
---|
589 | write(msg,'(a,i3)')'set_nfuel_cat: ifuelread=',ifuelread |
---|
590 | !$OMP END CRITICAL(SFIRE_PHYS_CRIT) |
---|
591 | call message(msg) |
---|
592 | |
---|
593 | if (ifuelread .eq. -1 .or. ifuelread .eq. 2) then |
---|
594 | !$OMP CRITICAL(SFIRE_PHYS_CRIT) |
---|
595 | call message('set_nfuel_cat: assuming nfuel_cat initialized already') |
---|
596 | call message(msg) |
---|
597 | !$OMP END CRITICAL(SFIRE_PHYS_CRIT) |
---|
598 | else if (ifuelread .eq. 0) then |
---|
599 | ! |
---|
600 | do j=jfts,jfte |
---|
601 | do i=ifts,ifte |
---|
602 | nfuel_cat(i,j)=real(nfuel_cat0) |
---|
603 | enddo |
---|
604 | enddo |
---|
605 | !$OMP CRITICAL(SFIRE_PHYS_CRIT) |
---|
606 | write(msg,'(a,i3)')'set_nfuel_cat: fuel initialized with category',nfuel_cat0 |
---|
607 | !$OMP END CRITICAL(SFIRE_PHYS_CRIT) |
---|
608 | call message(msg) |
---|
609 | |
---|
610 | else if (ifuelread .eq. 1) then |
---|
611 | ! |
---|
612 | ! make dependent on altitude (co mountains/forest vs. plains) |
---|
613 | ! 2000 m : 6562 ft ; 1600 m: 5249 ft |
---|
614 | |
---|
615 | ! ... user defines fuel category spatial variability ! param! |
---|
616 | do j=jfts,jfte |
---|
617 | do i=ifts,ifte |
---|
618 | ! nfuel_cat(i,j)= 2 ! grass with understory ! jm does nothing |
---|
619 | !jm t1=zsf(i,j)*slngth/100. |
---|
620 | t1 = zsf(i,j) ! this is in m |
---|
621 | if(t1.le.1524.)then ! up to 5000 ft |
---|
622 | nfuel_cat(i,j)= 3 ! tall grass |
---|
623 | else if(t1.ge.1524. .and. t1.le.2073.)then ! 5.0-6.8 kft. |
---|
624 | nfuel_cat(i,j)= 2 ! grass with understory |
---|
625 | else if(t1.ge.2073..and.t1.le.2438.)then ! 6.8-8.0 kft. |
---|
626 | nfuel_cat(i,j)= 8 ! timber litter - 10 (ponderosa) |
---|
627 | else if(t1.gt.2438. .and. t1.le. 3354.) then ! 8.0-11.0 kft. |
---|
628 | ! ... could also be mixed conifer. |
---|
629 | nfuel_cat(i,j)= 10 ! timber litter - 8 (lodgepole) |
---|
630 | else if(t1.gt.3354. .and. t1.le. 3658.) then ! 11.0-12.0 kft |
---|
631 | nfuel_cat(i,j)= 1 ! alpine meadow - 1 |
---|
632 | else if(t1.gt.3658. ) then ! > 12.0 kft |
---|
633 | nfuel_cat(i,j)= 14 ! no fuel. |
---|
634 | endif |
---|
635 | enddo |
---|
636 | enddo |
---|
637 | |
---|
638 | call message('set_nfuel_cat: fuel initialized by altitude') |
---|
639 | else |
---|
640 | |
---|
641 | call crash('set_nfuel_cat: bad ifuelread') |
---|
642 | endif |
---|
643 | ! .............end load fuel categories (or constant) here. |
---|
644 | |
---|
645 | end subroutine set_nfuel_cat |
---|
646 | |
---|
647 | ! |
---|
648 | !********************** |
---|
649 | ! |
---|
650 | |
---|
651 | subroutine fire_ros(ros_base,ros_wind,ros_slope, & |
---|
652 | propx,propy,i,j,fp) |
---|
653 | |
---|
654 | implicit none |
---|
655 | |
---|
656 | ! copied with the following changes ONLY: |
---|
657 | ! 0.5*(speed + abs(speed)) -> max(speed,0.) |
---|
658 | ! index l -> j |
---|
659 | ! not using nfuel_cat here - cell info was put into arrays passed as arguments |
---|
660 | ! in include file to avoid transcription errors when used elsewhere |
---|
661 | ! betaop is absorbed in phiwc, see module_fr_sfire_model/fire_startup |
---|
662 | ! return the base, wind, and slope contributions to the rate of spread separately |
---|
663 | ! because they may be needed to take advantage of known wind and slope vectors. |
---|
664 | ! They should add up to get the total rate of spread. |
---|
665 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
666 | ! ... calculates fire spread rate with McArthur formula or Rothermel |
---|
667 | ! using fuel type of fuel cell |
---|
668 | ! |
---|
669 | ! m/s =(ft/min) *.3048/60. =(ft/min) * .00508 ! conversion rate |
---|
670 | ! ft/min = m/s * 2.2369 * 88. = m/s * 196.850 ! conversion rate |
---|
671 | ! |
---|
672 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
673 | |
---|
674 | !*** arguments |
---|
675 | real, intent(out)::ros_base,ros_wind,ros_slope ! rate of spread contribution due to fuel, wind, and slope |
---|
676 | real, intent(in)::propx,propy |
---|
677 | integer, intent(in)::i,j ! node mesh coordinates |
---|
678 | type(fire_params),intent(in)::fp |
---|
679 | |
---|
680 | !*** local |
---|
681 | real:: speed, tanphi ! windspeed and slope in the direction normal to the fireline |
---|
682 | real:: umid, phis, phiw, spdms, umidm, excess |
---|
683 | real:: ros_back |
---|
684 | integer, parameter::ibeh=1 |
---|
685 | real, parameter::ros_max=6. |
---|
686 | character(len=128)msg |
---|
687 | real::cor_wind,cor_slope,nvx,nvy,scale |
---|
688 | |
---|
689 | |
---|
690 | !*** executable |
---|
691 | |
---|
692 | ! make sure normal direction is size 1 |
---|
693 | !scale=sqrt(propx*propx+propy*propy)+tiny(scale) |
---|
694 | scale=1. |
---|
695 | nvx=propx/scale |
---|
696 | nvy=propy/scale |
---|
697 | if (fire_advection.ne.0) then ! from flags in module_fr_sfire_util |
---|
698 | ! wind speed is total speed |
---|
699 | speed = sqrt(fp%vx(i,j)*fp%vx(i,j)+ fp%vy(i,j)*fp%vy(i,j))+tiny(speed) |
---|
700 | ! slope is total slope |
---|
701 | tanphi = sqrt(fp%dzdxf(i,j)*fp%dzdxf(i,j) + fp%dzdyf(i,j)*fp%dzdyf(i,j))+tiny(tanphi) |
---|
702 | ! cos of wind and spread, if >0 |
---|
703 | cor_wind = max(0.,(fp%vx(i,j)*nvx + fp%vy(i,j)*nvy)/speed) |
---|
704 | ! cos of slope and spread, if >0 |
---|
705 | cor_slope = max(0., (fp%dzdxf(i,j)*nvx + fp%dzdyf(i,j)*nvy)/tanphi) |
---|
706 | else |
---|
707 | ! wind speed in spread direction |
---|
708 | speed = fp%vx(i,j)*nvx + fp%vy(i,j)*nvy |
---|
709 | ! slope in spread direction |
---|
710 | tanphi = fp%dzdxf(i,j)*nvx + fp%dzdyf(i,j)*nvy |
---|
711 | cor_wind=1. |
---|
712 | cor_slope=1. |
---|
713 | endif |
---|
714 | |
---|
715 | if (.not. fp%ischap(i,j) > 0.) then ! if fuel is not chaparral, calculate rate of spread one of these ways |
---|
716 | if (ibeh .eq. 1) then ! use Rothermel formula |
---|
717 | ! ... if wind is 0 or into fireline, phiw = 0, &this reduces to backing ros. |
---|
718 | spdms = max(speed,0.) ! |
---|
719 | umidm = min(spdms,30.) ! max input wind spd is 30 m/s !param! |
---|
720 | umid = umidm * 196.850 ! m/s to ft/min |
---|
721 | ! eqn.: phiw = c * umid**bbb(i,j) * (fp%betafl(i,j)/betaop)**(-e) ! wind coef |
---|
722 | phiw = umid**fp%bbb(i,j) * fp%phiwc(i,j) ! wind coef |
---|
723 | phis=0. |
---|
724 | if (tanphi .gt. 0.) then |
---|
725 | phis = 5.275 *(fp%betafl(i,j))**(-0.3) *tanphi**2 ! slope factor |
---|
726 | endif |
---|
727 | ! rosm = fp%r_0(i,j)*(1. + phiw + phis) * .00508 ! spread rate, m/s |
---|
728 | ros_base = fp%r_0(i,j) * .00508 |
---|
729 | ros_wind = ros_base*phiw |
---|
730 | ros_slope= ros_base*phis |
---|
731 | |
---|
732 | |
---|
733 | else ! McArthur formula (Australian) |
---|
734 | ! rosm = 0.18*exp(0.8424*max(speed,0.)) |
---|
735 | ros_base = 0.18*exp(0.8424) |
---|
736 | ros_wind = 0.18*exp(0.8424*max(speed,0.)) |
---|
737 | ros_slope =0. |
---|
738 | endif |
---|
739 | ! |
---|
740 | else ! chaparral |
---|
741 | ! .... spread rate has no dependency on fuel character, only windspeed. |
---|
742 | spdms = max(speed,0.) |
---|
743 | ! rosm = 1.2974 * spdms**1.41 ! spread rate, m/s |
---|
744 | ! note: backing ros is 0 for chaparral without setting nozero value below |
---|
745 | !sp_n=.03333 |
---|
746 | ! chaparral backing fire spread rate 0.033 m/s ! param! |
---|
747 | !rosm= max(rosm, sp_n) ! no less than backing r.o.s. |
---|
748 | |
---|
749 | ros_back=.03333 ! chaparral backing fire spread rate 0.033 m/s ! param! |
---|
750 | ros_wind = 1.2974 * spdms**1.41 ! spread rate, m/s |
---|
751 | ros_wind = max(ros_wind, ros_back) |
---|
752 | ros_slope =0. |
---|
753 | |
---|
754 | endif |
---|
755 | ! if advection, multiply by the cosines |
---|
756 | ros_wind=ros_wind*cor_wind |
---|
757 | ros_slope=ros_slope*cor_slope |
---|
758 | ! |
---|
759 | ! ----------note! put an 6 m/s cap on max spread rate ----------- |
---|
760 | ! rosm= min(rosm, 6.) ! no faster than this cap ! param ! |
---|
761 | |
---|
762 | excess = ros_base + ros_wind + ros_slope - ros_max |
---|
763 | |
---|
764 | if (excess > 0.)then |
---|
765 | ! take it out of wind and slope in proportion |
---|
766 | ros_wind = ros_wind - excess*ros_wind/(ros_wind+ros_slope) |
---|
767 | ros_slope = ros_slope - excess*ros_slope/(ros_wind+ros_slope) |
---|
768 | endif |
---|
769 | |
---|
770 | !write(msg,*)i,j,' speed=',speed,' tanphi',tanphi,' ros=',ros_base,ros_wind,ros_slope |
---|
771 | !call message(msg) |
---|
772 | return |
---|
773 | |
---|
774 | contains |
---|
775 | real function nrm2(u,v) |
---|
776 | real, intent(in)::u,v |
---|
777 | nrm2=sqrt(u*u+v*v) |
---|
778 | end function nrm2 |
---|
779 | |
---|
780 | end subroutine fire_ros |
---|
781 | |
---|
782 | end module module_fr_sfire_phys |
---|