Skip to content

Spotting Models

These modules model firebrand generation and spot fire ignition. The Embers class implements the Albini (1979) physics-based model for firebrand lofting and ballistic transport from torching trees. The PerrymanSpotting class implements the Perryman et al. (2013) statistical model that samples firebrand landing locations from probability distributions parameterized by fireline intensity and wind speed.

The core simulation (FireSim) uses the physics-based Embers model, which tracks individual firebrand trajectories through the wind field. The prediction model (FirePredictor) uses the statistical PerrymanSpotting model, which is faster and better suited to the simplified prediction environment where fuel moisture is fixed and fire acceleration is not modeled.

Ember lofting, transport, and spot fire ignition model.

Simulate the generation of firebrands (embers) from torching trees, their ballistic transport through the atmosphere, and potential spot fire ignition upon landing. Based on the Albini (1979) torching-tree firebrand model.

Classes:

Name Description
- Embers

Manage firebrand lofting, flight, and landing for a simulation.

References

Albini, F. A. (1979). Spot fire distance from burning trees — a predictive model. USDA Forest Service General Technical Report INT-56.

Embers

Firebrand lofting, atmospheric transport, and spot fire ignition.

Manage a collection of airborne embers, lofting them from torching cells, transporting them through the wind field, and igniting spot fires where they land on fuel cells.

Attributes:

Name Type Description
ign_prob float

Probability of spot fire ignition per firebrand (0 to 1).

min_spot_dist_m float

Minimum horizontal distance (meters) a firebrand must travel before it can ignite a spot fire.

embers list[dict]

Currently airborne firebrands, each a dict with keys 'x', 'y', 'diam', 'height', 'start_elev', 'curr_time', 'elapsed'.

Source code in embrs/models/embers.py
 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
class Embers:
    """Firebrand lofting, atmospheric transport, and spot fire ignition.

    Manage a collection of airborne embers, lofting them from torching
    cells, transporting them through the wind field, and igniting spot
    fires where they land on fuel cells.

    Attributes:
        ign_prob (float): Probability of spot fire ignition per firebrand
            (0 to 1).
        min_spot_dist_m (float): Minimum horizontal distance (meters) a
            firebrand must travel before it can ignite a spot fire.
        embers (list[dict]): Currently airborne firebrands, each a dict
            with keys 'x', 'y', 'diam', 'height', 'start_elev',
            'curr_time', 'elapsed'.
    """

    def __init__(self, ign_prob: float, species: int, dbh: float,
                 min_spot_dist: float, limits: Tuple[float, float],
                 get_cell_from_xy: Callable[[float, float], Cell]):
        """Initialize the ember model.

        Args:
            ign_prob (float): Spot fire ignition probability per firebrand
                (0 to 1).
            species (int): Canopy species ID (see ``CanopySpecies``).
            dbh (float): Diameter at breast height of representative trees
                (centimeters).
            min_spot_dist (float): Minimum spotting distance (meters).
            limits (Tuple[float, float]): ``(x_lim, y_lim)`` simulation
                domain bounds (meters).
            get_cell_from_xy (Callable): Callback to retrieve a Cell from
                spatial coordinates.

        Raises:
            ValueError: If ``species`` is not in ``CanopySpecies``.
        """
        # Call back injection
        self.get_cell_from_xy = get_cell_from_xy

        self.x_lim, self.y_lim = limits

        # Spotting Ignition Probability
        self.ign_prob = ign_prob

        # Minimum spotting distance (prevents excessive spotting that is going to burn imminently anyways)
        self.min_spot_dist_m = min_spot_dist

        # Load the canopy species properties
        self.A = CanopySpecies.properties

        if CanopySpecies.species_names.get(species) is None:
            raise ValueError(f"Species {species} not found in CanopySpecies.")

        self.treespecies = species
        self.dbh = dbh # cm

        self.embers = []

    def loft(self, cell: Cell, sim_time_m: float):
        """Loft firebrands from a torching cell.

        Compute plume characteristics, then attempt to loft up to 16
        embers of increasing diameter. Each passes a probabilistic check
        against ``ign_prob``. Successfully lofted embers are appended to
        ``self.embers``.

        Args:
            cell (Cell): The burning cell that is torching.
            sim_time_m (float): Current simulation time (minutes).

        Side Effects:
            Sets ``cell.lofted = True``. Appends ember dicts to
            ``self.embers``.
        """
        cell.lofted = True

        self.plume(cell)
        crown_height_ft = m_to_ft(cell.canopy_height)
        z_0 = crown_height_ft

        # Attempts to loft 16 embers of different types
        for count in range(1, 17):
            diameter = 0.005 * count

            # Only loft embers that pass probability check
            prob_spot = np.random.random() # [0, 1]
            if prob_spot > self.ign_prob:
                continue

            zf = self.torchheight(diameter, z_0)
            if zf > 0:
                # Appends to self.embers
                ember = {
                    'x': cell.x_pos,
                    'y': cell.y_pos,
                    'source_x': cell.x_pos,
                    'source_y': cell.y_pos,
                    'diam': diameter,
                    'height': ft_to_m(zf),
                    'start_elev': cell.elevation_m,
                    'curr_time': sim_time_m, # Times are in minutes
                    'elapsed': 0.0
                }
                self.embers.append(ember)

            else:
                break

    def partcalc(self, vowf: float, z: float) -> float:
        """Compute time residual for iterating maximum particle height.

        Used by ``torchheight`` to solve for the maximum lofting height of
        a firebrand as a function of particle diameter and flame height
        (Albini 1979).

        Args:
            vowf (float): Normalized terminal velocity ratio (dimensionless).
            z (float): Candidate particle height (feet).

        Returns:
            float: Time residual (dimensionless). Converges to zero at the
                correct maximum height.
        """
        # Calculates tT for iterating max particle height z as a function of particle diameter and flameheight (steady_height)
        # for torching trees only, function used by torcheight
        a = 5.963
        b = 4.563
        r2 = (b + z / self.steady_height) / a
        r = np.sqrt(r2)
        tT1 = a / 3.0 * (r * r2 - 1.0)
        pt8vowf = 0.8 * vowf
        Temp = abs((1.0 - pt8vowf) / (1.0 - pt8vowf * r))
        tp = a / (pt8vowf**3) * (np.log(Temp) - pt8vowf * (r - 1.0) - 0.5 * pt8vowf**2 * (r2 - 1.0))
        return tp - tT1

    def torchheight(self, diameter: float, z_0: float) -> float:
        """Compute maximum lofting height for a firebrand (Albini 1979).

        Iteratively solve for the height a firebrand of given diameter
        can reach from a torching tree, accounting for plume buoyancy and
        particle terminal velocity.

        Args:
            diameter (float): Firebrand diameter (feet).
            z_0 (float): Crown height / initial release height (feet).

        Returns:
            float: Maximum lofting height (feet), or -1.0 if the firebrand
                cannot be lofted.
        """
        # Albini 1979 torching trees model
        vowf = 40.0 * np.sqrt(diameter / self.steady_height)
        if vowf < 1.0:
            zo_ratio = np.sqrt(z_0 / self.steady_height)
            if zo_ratio > vowf:
                Temp = abs((1.0 - vowf) / (zo_ratio - vowf))
                tf = 1.0 - zo_ratio + vowf * np.log(Temp)
                tt = 0.2 * vowf * (1.0 + vowf * np.log(1.0 + 1.0 / (1.0 - vowf)))
                aT = self.duration + 1.2 - tf - tt
                if aT >= 0.0:
                    inc = 2.0
                    inc2 = 1.0
                    z = inc * self.steady_height
                    tT = self.partcalc(vowf, z)
                    while abs(tT - aT) > 0.01:
                        if tT == aT:
                            break
                        else:
                            if tT < aT:
                                inc += inc2
                            else:
                                inc -= inc2
                                inc2 /= 2.0
                                inc += inc2
                        z = inc * self.steady_height
                        tT = self.partcalc(vowf, z)
                    return z
        return -1.0

    def calc_front_dist(self, cell: Cell) -> float:
        """Estimate the fire front length at a cell.

        Approximate the frontal distance as the cell side length times the
        number of burnable neighbors.

        Args:
            cell (Cell): Cell to evaluate.

        Returns:
            float: Estimated fire front length (meters).
        """
        a = cell.cell_size
        num_neighbors = len(cell.burnable_neighbors)

        dist = a * num_neighbors
        return dist

    def plume(self, cell: Cell):
        """Compute plume steady-state height and flame duration for a cell.

        Determine the plume characteristics based on crown fire status,
        crown fraction burned, canopy cover, and species properties. The
        number of torching trees is estimated from fire front length and
        canopy conditions.

        Args:
            cell (Cell): Cell providing crown status, CFB, and canopy data.

        Side Effects:
            Sets ``self.steady_height`` (feet) and ``self.duration``
            (dimensionless time).
        """
        treespecies = self.treespecies
        tnum = 1

        # DBH = diameter at breast height in cm
        # Convert DBH from cm to inches
        DBH = self.dbh / 2.54

        # compute frontal fire distance
        front_dist = self.calc_front_dist(cell)

        if cell._crown_status == CrownStatus.ACTIVE:
            # If the fire is active, use the active crown fire formula
            part1, part2 = self.A[treespecies][:2]
            tnum = max(1, int(front_dist / 5.0))
            self.steady_height = part1 * DBH**part2 * tnum**0.4
        else:
            if cell.cfb > 0.5 and cell.canopy_cover > 50:
                tnum += 1
            if cell.cfb > 0.8:
                tnum += 1
                if cell.canopy_cover > 50:
                    tnum = 6
                if cell.canopy_cover > 80:
                    tnum = 10

            part1, part2 = self.A[treespecies][:2]
            self.steady_height = part1 * DBH**part2 * tnum**0.4

        part3, part4 = self.A[treespecies][2:]
        self.duration = part3 * DBH**(-part4) * tnum**(-0.2)


    def vert_wind_speed(self, height_above_ground: float, canopy_ht: float,
                        wind_speed: float, wind_adj_factor: float) -> float:
        """Compute wind speed at a given height above ground.

        Below the canopy, return the midflame wind speed. Above the canopy,
        use the logarithmic wind profile (Albini & Baughman 1979).

        Args:
            height_above_ground (float): Height above ground (feet).
            canopy_ht (float): Canopy height (feet).
            wind_speed (float): 20-ft or reference wind speed (m/s).
            wind_adj_factor (float): Wind adjustment factor for sub-canopy.

        Returns:
            float: Wind speed at the specified height (ft/s).
        """
        if height_above_ground < canopy_ht:
            # Compute midflame wind speed
            midflame_wind  = wind_adj_factor * wind_speed
            u = m_to_ft(midflame_wind) # midflame wind speed in ft/s
        else:
            # Albini & Baughman 1979 Res. Pap. INT-221
            u = m_to_ft(wind_speed)/(np.log((20 + 0.36 * canopy_ht)/(0.1313 * canopy_ht)))

        return u

    def flight(self, end_curr_time_step: float) -> Set[Cell]:
        """Transport all airborne embers and ignite spot fires.

        Advance each ember through the wind field using a ballistic
        trajectory model until it either lands, exits the domain, or
        reaches the end of the current time step. Embers that land on
        fuel cells ignite spot fires.

        Args:
            end_curr_time_step (float): End time of the current simulation
                step (minutes).

        Returns:
            Set[Cell]: Set of cells where spot fires were ignited.

        Side Effects:
            Updates ``self.embers`` to retain only still-airborne embers.
            Sets ignited cells to ``CellStates.FIRE``.
        """
        spots = set()

        sstep = 0.25 # minutes

        carry = []
        for ember in self.embers:
            sx, sy = ember['x'], ember['y']

            diameter, Z, Zelev, curr_time, elapsed = ember['diam'], ember['height'], ember['start_elev'], ember['curr_time'], ember['elapsed']

            temp_ember = ember

            if curr_time == end_curr_time_step:
                carry.append(ember)
                continue
            else:
                curr_cell = self.get_cell_from_xy(sx, sy, oob_ok=True)                
                wind_speed, wind_dir = curr_cell.curr_wind()

                eH = m_to_ft(curr_cell.canopy_height)

                if eH <= 0.0:
                    eH = 1.0

                zfuel1 = curr_cell.elevation_m
                zfuel2 = curr_cell.elevation_m

                if elapsed == 0:
                    Zelev = zfuel1

                MAXZ = 39000*diameter

                voo = np.sqrt((1910.087*diameter)/0.18) # terminal velocity of particle in ft/s
                tao = (4.8 * voo) / (0.0064 * np.pi * 32)
                Xtot = 0 # reset total horizontal distance traveled
                Z = m_to_ft(Z) # convert to height in feet above original land
                eZelev = m_to_ft(Zelev)

                zt1 = 0.0
                rwinddir = wind_dir # degrees

                mZt = np.inf

                abort_flight = False
                while mZt > 1.0 and not abort_flight:
                    # Compute test vertical drop over ember over sstep
                    ztest1 = voo * tao * ((elapsed * 60)/tao - 0.5 * ((elapsed * 60)/tao) ** 2)
                    ztest2 = voo * tao * (((sstep + elapsed) * 60)/tao - 0.5 * (((sstep + elapsed) * 60)/tao) ** 2)

                    if ztest1 < ztest2: # Particle burned out before contacting ground
                        DZRate = (ztest2 - ztest1)/sstep

                    else:
                        tx = sx
                        ty = sy
                        break

                    passA = 0
                    while passA < 3:
                        last_cell = curr_cell

                        HtLimit = 50.0

                        passB = 0
                        while passB < 2:
                            subsstep = HtLimit/DZRate # Sub time step in minutes

                            if (curr_time + subsstep > end_curr_time_step):
                                subsstep = end_curr_time_step - curr_time
                                if subsstep < 0:
                                    subsstep = 0
                                    curr_time = end_curr_time_step
                                HtLimit = subsstep * DZRate # Recalc height limit in terms of new time step

                            zt2 = zt1 + HtLimit
                            if zt1 < 0:
                                Z1HtFromStart = eZelev-(Z-zt1) # Ember dropped below original land height
                            else:
                                Z1HtFromStart = eZelev+(Z-zt1)

                            if zt2 < 0:
                                Z2HtFromStart = eZelev-(Z-zt2) # Ember dropped below original land height
                            else:
                                Z2HtFromStart = eZelev+(Z-zt2)

                            Z1AboveGround = Z1HtFromStart - m_to_ft(zfuel1)
                            Z2AboveGround = Z2HtFromStart - m_to_ft(zfuel2)

                            mZt = ft_to_m(Z2AboveGround)

                            if mZt > 1.0:
                                wind_speed = wind_speed
                                wind_adj_factor = curr_cell.wind_adj_factor

                                UH = self.vert_wind_speed(Z1AboveGround, eH, wind_speed, wind_adj_factor)/2.0 # half of UeH in ft/s
                                UH += self.vert_wind_speed(Z2AboveGround, eH, wind_speed, wind_adj_factor)/2.0 # average UH in ft/s

                                if Z1AboveGround > 0:
                                    if Z2AboveGround > 0:
                                        ZAvgAboveGround = (Z1AboveGround + Z2AboveGround) / 2
                                    else:
                                        ZAvgAboveGround = Z1AboveGround / 2

                                    if ZAvgAboveGround > 0.1313 * eH:
                                        dxdt = UH/2.03*np.log(ZAvgAboveGround/(0.1313*eH)) # dxdt in ft/s
                                        Xt = dxdt * 60 * subsstep # ft (Albini 1979)
                                        mXt = ft_to_m(Xt) # incremental distance traveled m

                                    else:
                                        mXt = 0

                                    if passB < 1 and mXt > curr_cell.cell_size - 0.5: # Make sure we don't skip whole cells
                                        # Reduce the sub time step to avoid jumping a whole cell
                                        HtLimit /= 2
                                    else:
                                        break

                                passB += 1

                            else:
                                if mZt < -1.0:
                                    HtLimit = Z1AboveGround + (m_to_ft(zfuel1 - zfuel2))

                                    if HtLimit < 0: # Could happen if ember at zt1 is way below zfuel2 at zt2
                                        mXt = 0.0

                                        passA = 10 if elapsed == 0.0 else 6
                                        break

                                else:
                                    passA = 6
                                    break

                        if passA > 7:
                            break

                        # Update the position of the ember based on the wind direction
                        tx = sx + mXt * np.sin(np.deg2rad(rwinddir))
                        ty = sy + mXt * np.cos(np.deg2rad(rwinddir))

                        curr_cell = self.get_cell_from_xy(tx, ty, oob_ok = True)

                        # Make sure the ember doesn't land outside the sim region
                        if curr_cell is None:
                            # Ember is out of the map
                            passA = 10
                            abort_flight = True
                            break


                        if curr_cell != last_cell:
                            eH = m_to_ft(curr_cell.canopy_height)

                            if eH <= 0.0:
                                eH = 1.0

                            zfuel1 = zfuel2
                            zfuel2 = curr_cell.elevation_m
                            rwinddir = wind_dir

                        else:
                            if zt2 > MAXZ:
                                passA = 10
                                break
                            else:
                                passA +=2
                                break

                        passA += 1

                    if abort_flight:
                        break

                    sx = tx
                    sy = ty
                    elapsed += subsstep
                    Xtot += Xt
                    zt1 = zt2

                    if passA == 10:
                        # Cause failure of ember
                        curr_time = -1
                        mZt = 2.0
                        abort_flight = True
                        break
                    else:
                        if passA == 9:
                            mZt = 0.0
                        if curr_time < end_curr_time_step:
                            curr_time += subsstep
                        else:
                            break

                if np.abs(mZt) < 1.0 and curr_time < end_curr_time_step and not abort_flight:
                    if curr_cell.fuel.burnable:
                        prev_x = temp_ember['x']
                        prev_y = temp_ember['y']

                        dist = np.sqrt((prev_x - tx) ** 2 + (prev_y - ty) ** 2)

                        if dist > self.min_spot_dist_m:

                            if curr_cell.state == CellStates.FUEL:
                                curr_cell.spot_source_xy = (temp_ember['source_x'], temp_ember['source_y'])
                                spots.add(curr_cell) # Add to a set that will handle igniting spot fires
                                curr_cell.get_ign_params(0) # Get ignition parameters for the cell
                                curr_cell._set_state(CellStates.FIRE)
                else:
                    if curr_time == end_curr_time_step:
                        ember['x'] = sx
                        ember['y'] = sy
                        ember['curr_time'] = curr_time
                        ember['elapsed'] = elapsed

                        carry.append(ember)

        self.embers = carry

        return spots

__init__(ign_prob, species, dbh, min_spot_dist, limits, get_cell_from_xy)

Initialize the ember model.

Parameters:

Name Type Description Default
ign_prob float

Spot fire ignition probability per firebrand (0 to 1).

required
species int

Canopy species ID (see CanopySpecies).

required
dbh float

Diameter at breast height of representative trees (centimeters).

required
min_spot_dist float

Minimum spotting distance (meters).

required
limits Tuple[float, float]

(x_lim, y_lim) simulation domain bounds (meters).

required
get_cell_from_xy Callable

Callback to retrieve a Cell from spatial coordinates.

required

Raises:

Type Description
ValueError

If species is not in CanopySpecies.

Source code in embrs/models/embers.py
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
def __init__(self, ign_prob: float, species: int, dbh: float,
             min_spot_dist: float, limits: Tuple[float, float],
             get_cell_from_xy: Callable[[float, float], Cell]):
    """Initialize the ember model.

    Args:
        ign_prob (float): Spot fire ignition probability per firebrand
            (0 to 1).
        species (int): Canopy species ID (see ``CanopySpecies``).
        dbh (float): Diameter at breast height of representative trees
            (centimeters).
        min_spot_dist (float): Minimum spotting distance (meters).
        limits (Tuple[float, float]): ``(x_lim, y_lim)`` simulation
            domain bounds (meters).
        get_cell_from_xy (Callable): Callback to retrieve a Cell from
            spatial coordinates.

    Raises:
        ValueError: If ``species`` is not in ``CanopySpecies``.
    """
    # Call back injection
    self.get_cell_from_xy = get_cell_from_xy

    self.x_lim, self.y_lim = limits

    # Spotting Ignition Probability
    self.ign_prob = ign_prob

    # Minimum spotting distance (prevents excessive spotting that is going to burn imminently anyways)
    self.min_spot_dist_m = min_spot_dist

    # Load the canopy species properties
    self.A = CanopySpecies.properties

    if CanopySpecies.species_names.get(species) is None:
        raise ValueError(f"Species {species} not found in CanopySpecies.")

    self.treespecies = species
    self.dbh = dbh # cm

    self.embers = []

calc_front_dist(cell)

Estimate the fire front length at a cell.

Approximate the frontal distance as the cell side length times the number of burnable neighbors.

Parameters:

Name Type Description Default
cell Cell

Cell to evaluate.

required

Returns:

Name Type Description
float float

Estimated fire front length (meters).

Source code in embrs/models/embers.py
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
def calc_front_dist(self, cell: Cell) -> float:
    """Estimate the fire front length at a cell.

    Approximate the frontal distance as the cell side length times the
    number of burnable neighbors.

    Args:
        cell (Cell): Cell to evaluate.

    Returns:
        float: Estimated fire front length (meters).
    """
    a = cell.cell_size
    num_neighbors = len(cell.burnable_neighbors)

    dist = a * num_neighbors
    return dist

flight(end_curr_time_step)

Transport all airborne embers and ignite spot fires.

Advance each ember through the wind field using a ballistic trajectory model until it either lands, exits the domain, or reaches the end of the current time step. Embers that land on fuel cells ignite spot fires.

Parameters:

Name Type Description Default
end_curr_time_step float

End time of the current simulation step (minutes).

required

Returns:

Type Description
Set[Cell]

Set[Cell]: Set of cells where spot fires were ignited.

Side Effects

Updates self.embers to retain only still-airborne embers. Sets ignited cells to CellStates.FIRE.

Source code in embrs/models/embers.py
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
def flight(self, end_curr_time_step: float) -> Set[Cell]:
    """Transport all airborne embers and ignite spot fires.

    Advance each ember through the wind field using a ballistic
    trajectory model until it either lands, exits the domain, or
    reaches the end of the current time step. Embers that land on
    fuel cells ignite spot fires.

    Args:
        end_curr_time_step (float): End time of the current simulation
            step (minutes).

    Returns:
        Set[Cell]: Set of cells where spot fires were ignited.

    Side Effects:
        Updates ``self.embers`` to retain only still-airborne embers.
        Sets ignited cells to ``CellStates.FIRE``.
    """
    spots = set()

    sstep = 0.25 # minutes

    carry = []
    for ember in self.embers:
        sx, sy = ember['x'], ember['y']

        diameter, Z, Zelev, curr_time, elapsed = ember['diam'], ember['height'], ember['start_elev'], ember['curr_time'], ember['elapsed']

        temp_ember = ember

        if curr_time == end_curr_time_step:
            carry.append(ember)
            continue
        else:
            curr_cell = self.get_cell_from_xy(sx, sy, oob_ok=True)                
            wind_speed, wind_dir = curr_cell.curr_wind()

            eH = m_to_ft(curr_cell.canopy_height)

            if eH <= 0.0:
                eH = 1.0

            zfuel1 = curr_cell.elevation_m
            zfuel2 = curr_cell.elevation_m

            if elapsed == 0:
                Zelev = zfuel1

            MAXZ = 39000*diameter

            voo = np.sqrt((1910.087*diameter)/0.18) # terminal velocity of particle in ft/s
            tao = (4.8 * voo) / (0.0064 * np.pi * 32)
            Xtot = 0 # reset total horizontal distance traveled
            Z = m_to_ft(Z) # convert to height in feet above original land
            eZelev = m_to_ft(Zelev)

            zt1 = 0.0
            rwinddir = wind_dir # degrees

            mZt = np.inf

            abort_flight = False
            while mZt > 1.0 and not abort_flight:
                # Compute test vertical drop over ember over sstep
                ztest1 = voo * tao * ((elapsed * 60)/tao - 0.5 * ((elapsed * 60)/tao) ** 2)
                ztest2 = voo * tao * (((sstep + elapsed) * 60)/tao - 0.5 * (((sstep + elapsed) * 60)/tao) ** 2)

                if ztest1 < ztest2: # Particle burned out before contacting ground
                    DZRate = (ztest2 - ztest1)/sstep

                else:
                    tx = sx
                    ty = sy
                    break

                passA = 0
                while passA < 3:
                    last_cell = curr_cell

                    HtLimit = 50.0

                    passB = 0
                    while passB < 2:
                        subsstep = HtLimit/DZRate # Sub time step in minutes

                        if (curr_time + subsstep > end_curr_time_step):
                            subsstep = end_curr_time_step - curr_time
                            if subsstep < 0:
                                subsstep = 0
                                curr_time = end_curr_time_step
                            HtLimit = subsstep * DZRate # Recalc height limit in terms of new time step

                        zt2 = zt1 + HtLimit
                        if zt1 < 0:
                            Z1HtFromStart = eZelev-(Z-zt1) # Ember dropped below original land height
                        else:
                            Z1HtFromStart = eZelev+(Z-zt1)

                        if zt2 < 0:
                            Z2HtFromStart = eZelev-(Z-zt2) # Ember dropped below original land height
                        else:
                            Z2HtFromStart = eZelev+(Z-zt2)

                        Z1AboveGround = Z1HtFromStart - m_to_ft(zfuel1)
                        Z2AboveGround = Z2HtFromStart - m_to_ft(zfuel2)

                        mZt = ft_to_m(Z2AboveGround)

                        if mZt > 1.0:
                            wind_speed = wind_speed
                            wind_adj_factor = curr_cell.wind_adj_factor

                            UH = self.vert_wind_speed(Z1AboveGround, eH, wind_speed, wind_adj_factor)/2.0 # half of UeH in ft/s
                            UH += self.vert_wind_speed(Z2AboveGround, eH, wind_speed, wind_adj_factor)/2.0 # average UH in ft/s

                            if Z1AboveGround > 0:
                                if Z2AboveGround > 0:
                                    ZAvgAboveGround = (Z1AboveGround + Z2AboveGround) / 2
                                else:
                                    ZAvgAboveGround = Z1AboveGround / 2

                                if ZAvgAboveGround > 0.1313 * eH:
                                    dxdt = UH/2.03*np.log(ZAvgAboveGround/(0.1313*eH)) # dxdt in ft/s
                                    Xt = dxdt * 60 * subsstep # ft (Albini 1979)
                                    mXt = ft_to_m(Xt) # incremental distance traveled m

                                else:
                                    mXt = 0

                                if passB < 1 and mXt > curr_cell.cell_size - 0.5: # Make sure we don't skip whole cells
                                    # Reduce the sub time step to avoid jumping a whole cell
                                    HtLimit /= 2
                                else:
                                    break

                            passB += 1

                        else:
                            if mZt < -1.0:
                                HtLimit = Z1AboveGround + (m_to_ft(zfuel1 - zfuel2))

                                if HtLimit < 0: # Could happen if ember at zt1 is way below zfuel2 at zt2
                                    mXt = 0.0

                                    passA = 10 if elapsed == 0.0 else 6
                                    break

                            else:
                                passA = 6
                                break

                    if passA > 7:
                        break

                    # Update the position of the ember based on the wind direction
                    tx = sx + mXt * np.sin(np.deg2rad(rwinddir))
                    ty = sy + mXt * np.cos(np.deg2rad(rwinddir))

                    curr_cell = self.get_cell_from_xy(tx, ty, oob_ok = True)

                    # Make sure the ember doesn't land outside the sim region
                    if curr_cell is None:
                        # Ember is out of the map
                        passA = 10
                        abort_flight = True
                        break


                    if curr_cell != last_cell:
                        eH = m_to_ft(curr_cell.canopy_height)

                        if eH <= 0.0:
                            eH = 1.0

                        zfuel1 = zfuel2
                        zfuel2 = curr_cell.elevation_m
                        rwinddir = wind_dir

                    else:
                        if zt2 > MAXZ:
                            passA = 10
                            break
                        else:
                            passA +=2
                            break

                    passA += 1

                if abort_flight:
                    break

                sx = tx
                sy = ty
                elapsed += subsstep
                Xtot += Xt
                zt1 = zt2

                if passA == 10:
                    # Cause failure of ember
                    curr_time = -1
                    mZt = 2.0
                    abort_flight = True
                    break
                else:
                    if passA == 9:
                        mZt = 0.0
                    if curr_time < end_curr_time_step:
                        curr_time += subsstep
                    else:
                        break

            if np.abs(mZt) < 1.0 and curr_time < end_curr_time_step and not abort_flight:
                if curr_cell.fuel.burnable:
                    prev_x = temp_ember['x']
                    prev_y = temp_ember['y']

                    dist = np.sqrt((prev_x - tx) ** 2 + (prev_y - ty) ** 2)

                    if dist > self.min_spot_dist_m:

                        if curr_cell.state == CellStates.FUEL:
                            curr_cell.spot_source_xy = (temp_ember['source_x'], temp_ember['source_y'])
                            spots.add(curr_cell) # Add to a set that will handle igniting spot fires
                            curr_cell.get_ign_params(0) # Get ignition parameters for the cell
                            curr_cell._set_state(CellStates.FIRE)
            else:
                if curr_time == end_curr_time_step:
                    ember['x'] = sx
                    ember['y'] = sy
                    ember['curr_time'] = curr_time
                    ember['elapsed'] = elapsed

                    carry.append(ember)

    self.embers = carry

    return spots

loft(cell, sim_time_m)

Loft firebrands from a torching cell.

Compute plume characteristics, then attempt to loft up to 16 embers of increasing diameter. Each passes a probabilistic check against ign_prob. Successfully lofted embers are appended to self.embers.

Parameters:

Name Type Description Default
cell Cell

The burning cell that is torching.

required
sim_time_m float

Current simulation time (minutes).

required
Side Effects

Sets cell.lofted = True. Appends ember dicts to self.embers.

Source code in embrs/models/embers.py
 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
def loft(self, cell: Cell, sim_time_m: float):
    """Loft firebrands from a torching cell.

    Compute plume characteristics, then attempt to loft up to 16
    embers of increasing diameter. Each passes a probabilistic check
    against ``ign_prob``. Successfully lofted embers are appended to
    ``self.embers``.

    Args:
        cell (Cell): The burning cell that is torching.
        sim_time_m (float): Current simulation time (minutes).

    Side Effects:
        Sets ``cell.lofted = True``. Appends ember dicts to
        ``self.embers``.
    """
    cell.lofted = True

    self.plume(cell)
    crown_height_ft = m_to_ft(cell.canopy_height)
    z_0 = crown_height_ft

    # Attempts to loft 16 embers of different types
    for count in range(1, 17):
        diameter = 0.005 * count

        # Only loft embers that pass probability check
        prob_spot = np.random.random() # [0, 1]
        if prob_spot > self.ign_prob:
            continue

        zf = self.torchheight(diameter, z_0)
        if zf > 0:
            # Appends to self.embers
            ember = {
                'x': cell.x_pos,
                'y': cell.y_pos,
                'source_x': cell.x_pos,
                'source_y': cell.y_pos,
                'diam': diameter,
                'height': ft_to_m(zf),
                'start_elev': cell.elevation_m,
                'curr_time': sim_time_m, # Times are in minutes
                'elapsed': 0.0
            }
            self.embers.append(ember)

        else:
            break

partcalc(vowf, z)

Compute time residual for iterating maximum particle height.

Used by torchheight to solve for the maximum lofting height of a firebrand as a function of particle diameter and flame height (Albini 1979).

Parameters:

Name Type Description Default
vowf float

Normalized terminal velocity ratio (dimensionless).

required
z float

Candidate particle height (feet).

required

Returns:

Name Type Description
float float

Time residual (dimensionless). Converges to zero at the correct maximum height.

Source code in embrs/models/embers.py
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
def partcalc(self, vowf: float, z: float) -> float:
    """Compute time residual for iterating maximum particle height.

    Used by ``torchheight`` to solve for the maximum lofting height of
    a firebrand as a function of particle diameter and flame height
    (Albini 1979).

    Args:
        vowf (float): Normalized terminal velocity ratio (dimensionless).
        z (float): Candidate particle height (feet).

    Returns:
        float: Time residual (dimensionless). Converges to zero at the
            correct maximum height.
    """
    # Calculates tT for iterating max particle height z as a function of particle diameter and flameheight (steady_height)
    # for torching trees only, function used by torcheight
    a = 5.963
    b = 4.563
    r2 = (b + z / self.steady_height) / a
    r = np.sqrt(r2)
    tT1 = a / 3.0 * (r * r2 - 1.0)
    pt8vowf = 0.8 * vowf
    Temp = abs((1.0 - pt8vowf) / (1.0 - pt8vowf * r))
    tp = a / (pt8vowf**3) * (np.log(Temp) - pt8vowf * (r - 1.0) - 0.5 * pt8vowf**2 * (r2 - 1.0))
    return tp - tT1

plume(cell)

Compute plume steady-state height and flame duration for a cell.

Determine the plume characteristics based on crown fire status, crown fraction burned, canopy cover, and species properties. The number of torching trees is estimated from fire front length and canopy conditions.

Parameters:

Name Type Description Default
cell Cell

Cell providing crown status, CFB, and canopy data.

required
Side Effects

Sets self.steady_height (feet) and self.duration (dimensionless time).

Source code in embrs/models/embers.py
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
def plume(self, cell: Cell):
    """Compute plume steady-state height and flame duration for a cell.

    Determine the plume characteristics based on crown fire status,
    crown fraction burned, canopy cover, and species properties. The
    number of torching trees is estimated from fire front length and
    canopy conditions.

    Args:
        cell (Cell): Cell providing crown status, CFB, and canopy data.

    Side Effects:
        Sets ``self.steady_height`` (feet) and ``self.duration``
        (dimensionless time).
    """
    treespecies = self.treespecies
    tnum = 1

    # DBH = diameter at breast height in cm
    # Convert DBH from cm to inches
    DBH = self.dbh / 2.54

    # compute frontal fire distance
    front_dist = self.calc_front_dist(cell)

    if cell._crown_status == CrownStatus.ACTIVE:
        # If the fire is active, use the active crown fire formula
        part1, part2 = self.A[treespecies][:2]
        tnum = max(1, int(front_dist / 5.0))
        self.steady_height = part1 * DBH**part2 * tnum**0.4
    else:
        if cell.cfb > 0.5 and cell.canopy_cover > 50:
            tnum += 1
        if cell.cfb > 0.8:
            tnum += 1
            if cell.canopy_cover > 50:
                tnum = 6
            if cell.canopy_cover > 80:
                tnum = 10

        part1, part2 = self.A[treespecies][:2]
        self.steady_height = part1 * DBH**part2 * tnum**0.4

    part3, part4 = self.A[treespecies][2:]
    self.duration = part3 * DBH**(-part4) * tnum**(-0.2)

torchheight(diameter, z_0)

Compute maximum lofting height for a firebrand (Albini 1979).

Iteratively solve for the height a firebrand of given diameter can reach from a torching tree, accounting for plume buoyancy and particle terminal velocity.

Parameters:

Name Type Description Default
diameter float

Firebrand diameter (feet).

required
z_0 float

Crown height / initial release height (feet).

required

Returns:

Name Type Description
float float

Maximum lofting height (feet), or -1.0 if the firebrand cannot be lofted.

Source code in embrs/models/embers.py
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
def torchheight(self, diameter: float, z_0: float) -> float:
    """Compute maximum lofting height for a firebrand (Albini 1979).

    Iteratively solve for the height a firebrand of given diameter
    can reach from a torching tree, accounting for plume buoyancy and
    particle terminal velocity.

    Args:
        diameter (float): Firebrand diameter (feet).
        z_0 (float): Crown height / initial release height (feet).

    Returns:
        float: Maximum lofting height (feet), or -1.0 if the firebrand
            cannot be lofted.
    """
    # Albini 1979 torching trees model
    vowf = 40.0 * np.sqrt(diameter / self.steady_height)
    if vowf < 1.0:
        zo_ratio = np.sqrt(z_0 / self.steady_height)
        if zo_ratio > vowf:
            Temp = abs((1.0 - vowf) / (zo_ratio - vowf))
            tf = 1.0 - zo_ratio + vowf * np.log(Temp)
            tt = 0.2 * vowf * (1.0 + vowf * np.log(1.0 + 1.0 / (1.0 - vowf)))
            aT = self.duration + 1.2 - tf - tt
            if aT >= 0.0:
                inc = 2.0
                inc2 = 1.0
                z = inc * self.steady_height
                tT = self.partcalc(vowf, z)
                while abs(tT - aT) > 0.01:
                    if tT == aT:
                        break
                    else:
                        if tT < aT:
                            inc += inc2
                        else:
                            inc -= inc2
                            inc2 /= 2.0
                            inc += inc2
                    z = inc * self.steady_height
                    tT = self.partcalc(vowf, z)
                return z
    return -1.0

vert_wind_speed(height_above_ground, canopy_ht, wind_speed, wind_adj_factor)

Compute wind speed at a given height above ground.

Below the canopy, return the midflame wind speed. Above the canopy, use the logarithmic wind profile (Albini & Baughman 1979).

Parameters:

Name Type Description Default
height_above_ground float

Height above ground (feet).

required
canopy_ht float

Canopy height (feet).

required
wind_speed float

20-ft or reference wind speed (m/s).

required
wind_adj_factor float

Wind adjustment factor for sub-canopy.

required

Returns:

Name Type Description
float float

Wind speed at the specified height (ft/s).

Source code in embrs/models/embers.py
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
def vert_wind_speed(self, height_above_ground: float, canopy_ht: float,
                    wind_speed: float, wind_adj_factor: float) -> float:
    """Compute wind speed at a given height above ground.

    Below the canopy, return the midflame wind speed. Above the canopy,
    use the logarithmic wind profile (Albini & Baughman 1979).

    Args:
        height_above_ground (float): Height above ground (feet).
        canopy_ht (float): Canopy height (feet).
        wind_speed (float): 20-ft or reference wind speed (m/s).
        wind_adj_factor (float): Wind adjustment factor for sub-canopy.

    Returns:
        float: Wind speed at the specified height (ft/s).
    """
    if height_above_ground < canopy_ht:
        # Compute midflame wind speed
        midflame_wind  = wind_adj_factor * wind_speed
        u = m_to_ft(midflame_wind) # midflame wind speed in ft/s
    else:
        # Albini & Baughman 1979 Res. Pap. INT-221
        u = m_to_ft(wind_speed)/(np.log((20 + 0.36 * canopy_ht)/(0.1313 * canopy_ht)))

    return u

BTU_ft2_min_to_kW_m2(f_btu_ft2_min)

Convert heat flux from BTU/(ft^2*min) to kW/m^2.

Parameters:

Name Type Description Default
f_btu_ft2_min float

Heat flux in BTU/(ft^2*min).

required

Returns:

Name Type Description
float float

Heat flux in kW/m^2.

Source code in embrs/utilities/unit_conversions.py
220
221
222
223
224
225
226
227
228
229
def BTU_ft2_min_to_kW_m2(f_btu_ft2_min: float) -> float:
    """Convert heat flux from BTU/(ft^2*min) to kW/m^2.

    Args:
        f_btu_ft2_min (float): Heat flux in BTU/(ft^2*min).

    Returns:
        float: Heat flux in kW/m^2.
    """
    return f_btu_ft2_min * _BTU_FT2_MIN_TO_KW_M2

BTU_ft_min_to_kW_m(f_btu_ft_min)

Convert fireline intensity from BTU/(ft*min) to kW/m.

Parameters:

Name Type Description Default
f_btu_ft_min float

Fireline intensity in BTU/(ft*min).

required

Returns:

Name Type Description
float float

Fireline intensity in kW/m.

Source code in embrs/utilities/unit_conversions.py
232
233
234
235
236
237
238
239
240
241
def BTU_ft_min_to_kW_m(f_btu_ft_min: float) -> float:
    """Convert fireline intensity from BTU/(ft*min) to kW/m.

    Args:
        f_btu_ft_min (float): Fireline intensity in BTU/(ft*min).

    Returns:
        float: Fireline intensity in kW/m.
    """
    return f_btu_ft_min * _BTU_FT_MIN_TO_KW_M

BTU_ft_min_to_kcal_s_m(f_btu_ft_min)

Convert fireline intensity from BTU/(ftmin) to kcal/(sm).

Parameters:

Name Type Description Default
f_btu_ft_min float

Fireline intensity in BTU/(ft*min).

required

Returns:

Name Type Description
float float

Fireline intensity in kcal/(s*m).

Source code in embrs/utilities/unit_conversions.py
244
245
246
247
248
249
250
251
252
253
def BTU_ft_min_to_kcal_s_m(f_btu_ft_min: float) -> float:
    """Convert fireline intensity from BTU/(ft*min) to kcal/(s*m).

    Args:
        f_btu_ft_min (float): Fireline intensity in BTU/(ft*min).

    Returns:
        float: Fireline intensity in kcal/(s*m).
    """
    return f_btu_ft_min * _BTU_FT_MIN_TO_KCAL_S_M

BTU_lb_to_cal_g(f_btu_lb)

Convert heat content from BTU/lb to cal/g.

Parameters:

Name Type Description Default
f_btu_lb float

Heat content in BTU/lb.

required

Returns:

Name Type Description
float float

Heat content in cal/g.

Source code in embrs/utilities/unit_conversions.py
272
273
274
275
276
277
278
279
280
281
def BTU_lb_to_cal_g(f_btu_lb: float) -> float:
    """Convert heat content from BTU/lb to cal/g.

    Args:
        f_btu_lb (float): Heat content in BTU/lb.

    Returns:
        float: Heat content in cal/g.
    """
    return f_btu_lb * _BTU_LB_TO_CAL_G

F_to_C(f_f)

Convert temperature from Fahrenheit to Celsius.

Parameters:

Name Type Description Default
f_f float

Temperature in degrees Fahrenheit.

required

Returns:

Name Type Description
float float

Temperature in degrees Celsius.

Source code in embrs/utilities/unit_conversions.py
48
49
50
51
52
53
54
55
56
57
def F_to_C(f_f: float) -> float:
    """Convert temperature from Fahrenheit to Celsius.

    Args:
        f_f (float): Temperature in degrees Fahrenheit.

    Returns:
        float: Temperature in degrees Celsius.
    """
    return (5.0 / 9.0) * (f_f - 32.0)

KiSq_to_Lbsft2(f_kisq)

Convert fuel loading from kg/m^2 to lb/ft^2.

Parameters:

Name Type Description Default
f_kisq float

Fuel loading in kg/m^2.

required

Returns:

Name Type Description
float float

Fuel loading in lb/ft^2.

Source code in embrs/utilities/unit_conversions.py
156
157
158
159
160
161
162
163
164
165
def KiSq_to_Lbsft2(f_kisq: float) -> float:
    """Convert fuel loading from kg/m^2 to lb/ft^2.

    Args:
        f_kisq (float): Fuel loading in kg/m^2.

    Returns:
        float: Fuel loading in lb/ft^2.
    """
    return f_kisq * _KISQ_TO_LBSFT2

KiSq_to_TPA(f_kisq)

Convert fuel loading from kg/m^2 to tons per acre.

Parameters:

Name Type Description Default
f_kisq float

Fuel loading in kg/m^2.

required

Returns:

Name Type Description
float float

Fuel loading in tons per acre.

Source code in embrs/utilities/unit_conversions.py
204
205
206
207
208
209
210
211
212
213
def KiSq_to_TPA(f_kisq: float) -> float:
    """Convert fuel loading from kg/m^2 to tons per acre.

    Args:
        f_kisq (float): Fuel loading in kg/m^2.

    Returns:
        float: Fuel loading in tons per acre.
    """
    return f_kisq * _KISQ_TO_TPA

Lbsft2_to_KiSq(f_libsft2)

Convert fuel loading from lb/ft^2 to kg/m^2.

Parameters:

Name Type Description Default
f_libsft2 float

Fuel loading in lb/ft^2.

required

Returns:

Name Type Description
float float

Fuel loading in kg/m^2.

Source code in embrs/utilities/unit_conversions.py
144
145
146
147
148
149
150
151
152
153
def Lbsft2_to_KiSq(f_libsft2: float) -> float:
    """Convert fuel loading from lb/ft^2 to kg/m^2.

    Args:
        f_libsft2 (float): Fuel loading in lb/ft^2.

    Returns:
        float: Fuel loading in kg/m^2.
    """
    return f_libsft2 * _LBSFT2_TO_KISQ

Lbsft2_to_TPA(f_lbsft2)

Convert fuel loading from lb/ft^2 to tons per acre.

Parameters:

Name Type Description Default
f_lbsft2 float

Fuel loading in lb/ft^2.

required

Returns:

Name Type Description
float float

Fuel loading in tons per acre.

Source code in embrs/utilities/unit_conversions.py
192
193
194
195
196
197
198
199
200
201
def Lbsft2_to_TPA(f_lbsft2: float) -> float:
    """Convert fuel loading from lb/ft^2 to tons per acre.

    Args:
        f_lbsft2 (float): Fuel loading in lb/ft^2.

    Returns:
        float: Fuel loading in tons per acre.
    """
    return f_lbsft2 * _LBSFT2_TO_TPA

TPA_to_KiSq(f_tpa)

Convert fuel loading from tons per acre to kg/m^2.

Parameters:

Name Type Description Default
f_tpa float

Fuel loading in tons per acre.

required

Returns:

Name Type Description
float float

Fuel loading in kg/m^2.

Source code in embrs/utilities/unit_conversions.py
168
169
170
171
172
173
174
175
176
177
def TPA_to_KiSq(f_tpa: float) -> float:
    """Convert fuel loading from tons per acre to kg/m^2.

    Args:
        f_tpa (float): Fuel loading in tons per acre.

    Returns:
        float: Fuel loading in kg/m^2.
    """
    return f_tpa / _TPA_TO_KISQ_DIVISOR

TPA_to_Lbsft2(f_tpa)

Convert fuel loading from tons per acre to lb/ft^2.

Parameters:

Name Type Description Default
f_tpa float

Fuel loading in tons per acre.

required

Returns:

Name Type Description
float float

Fuel loading in lb/ft^2.

Source code in embrs/utilities/unit_conversions.py
180
181
182
183
184
185
186
187
188
189
def TPA_to_Lbsft2(f_tpa: float) -> float:
    """Convert fuel loading from tons per acre to lb/ft^2.

    Args:
        f_tpa (float): Fuel loading in tons per acre.

    Returns:
        float: Fuel loading in lb/ft^2.
    """
    return f_tpa * _TPA_TO_LBSFT2

cal_g_to_BTU_lb(f_cal_g)

Convert heat content from cal/g to BTU/lb.

Parameters:

Name Type Description Default
f_cal_g float

Heat content in cal/g.

required

Returns:

Name Type Description
float float

Heat content in BTU/lb.

Source code in embrs/utilities/unit_conversions.py
260
261
262
263
264
265
266
267
268
269
def cal_g_to_BTU_lb(f_cal_g: float) -> float:
    """Convert heat content from cal/g to BTU/lb.

    Args:
        f_cal_g (float): Heat content in cal/g.

    Returns:
        float: Heat content in BTU/lb.
    """
    return f_cal_g * _CAL_G_TO_BTU_LB

ft_min_to_m_s(f_ft_min)

Convert speed from feet per minute to meters per second.

Parameters:

Name Type Description Default
f_ft_min float

Speed in ft/min.

required

Returns:

Name Type Description
float float

Speed in m/s.

Source code in embrs/utilities/unit_conversions.py
 92
 93
 94
 95
 96
 97
 98
 99
100
101
def ft_min_to_m_s(f_ft_min: float) -> float:
    """Convert speed from feet per minute to meters per second.

    Args:
        f_ft_min (float): Speed in ft/min.

    Returns:
        float: Speed in m/s.
    """
    return f_ft_min * _FT_MIN_TO_M_S

ft_min_to_mph(f_ft_min)

Convert speed from feet per minute to miles per hour.

Parameters:

Name Type Description Default
f_ft_min float

Speed in ft/min.

required

Returns:

Name Type Description
float float

Speed in mph.

Source code in embrs/utilities/unit_conversions.py
116
117
118
119
120
121
122
123
124
125
def ft_min_to_mph(f_ft_min: float) -> float:
    """Convert speed from feet per minute to miles per hour.

    Args:
        f_ft_min (float): Speed in ft/min.

    Returns:
        float: Speed in mph.
    """
    return f_ft_min * _FT_MIN_TO_MPH

ft_to_m(f_ft)

Convert length from feet to meters.

Parameters:

Name Type Description Default
f_ft float

Length in feet.

required

Returns:

Name Type Description
float float

Length in meters.

Source code in embrs/utilities/unit_conversions.py
76
77
78
79
80
81
82
83
84
85
def ft_to_m(f_ft: float) -> float:
    """Convert length from feet to meters.

    Args:
        f_ft (float): Length in feet.

    Returns:
        float: Length in meters.
    """
    return f_ft * _FT_TO_M

m_s_to_ft_min(m_s)

Convert speed from meters per second to feet per minute.

Parameters:

Name Type Description Default
m_s float

Speed in m/s.

required

Returns:

Name Type Description
float float

Speed in ft/min.

Source code in embrs/utilities/unit_conversions.py
104
105
106
107
108
109
110
111
112
113
def m_s_to_ft_min(m_s: float) -> float:
    """Convert speed from meters per second to feet per minute.

    Args:
        m_s (float): Speed in m/s.

    Returns:
        float: Speed in ft/min.
    """
    return m_s * _M_S_TO_FT_MIN

m_to_ft(f_m)

Convert length from meters to feet.

Parameters:

Name Type Description Default
f_m float

Length in meters.

required

Returns:

Name Type Description
float float

Length in feet.

Source code in embrs/utilities/unit_conversions.py
64
65
66
67
68
69
70
71
72
73
def m_to_ft(f_m: float) -> float:
    """Convert length from meters to feet.

    Args:
        f_m (float): Length in meters.

    Returns:
        float: Length in feet.
    """
    return f_m * _M_TO_FT

mph_to_ft_min(f_mph)

Convert speed from miles per hour to feet per minute.

Parameters:

Name Type Description Default
f_mph float

Speed in mph.

required

Returns:

Name Type Description
float float

Speed in ft/min.

Source code in embrs/utilities/unit_conversions.py
128
129
130
131
132
133
134
135
136
137
def mph_to_ft_min(f_mph: float) -> float:
    """Convert speed from miles per hour to feet per minute.

    Args:
        f_mph (float): Speed in mph.

    Returns:
        float: Speed in ft/min.
    """
    return f_mph * _MPH_TO_FT_MIN

Statistical firebrand spotting model based on Perryman et al.

Generate firebrand landing locations using lognormal and normal distributions parameterized by fireline intensity and wind speed. The parallel (downwind) distribution depends on the Froude number regime; the perpendicular distribution is scaled by cell width.

Classes:

Name Description
- PerrymanSpotting

Sample firebrand landing positions from a burning cell.

References

Perryman, H. A., et al. (2013). A mathematical model of spot fire transport. International Journal of Wildland Fire, 22, 350-358.

PerrymanSpotting

Statistical firebrand landing model.

Sample firebrand landing coordinates from probability distributions parameterized by wind speed and fireline intensity. Firebrands that land within 50 meters of the source cell are filtered out.

Attributes:

Name Type Description
embers list[dict]

Accumulated firebrand landing positions, each with keys 'x', 'y', 'd' (distance in meters).

spot_delay_s float

Delay before spotted fires can ignite (seconds).

Source code in embrs/models/perryman_spot.py
 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
class PerrymanSpotting:
    """Statistical firebrand landing model.

    Sample firebrand landing coordinates from probability distributions
    parameterized by wind speed and fireline intensity. Firebrands that
    land within 50 meters of the source cell are filtered out.

    Attributes:
        embers (list[dict]): Accumulated firebrand landing positions,
            each with keys 'x', 'y', 'd' (distance in meters).
        spot_delay_s (float): Delay before spotted fires can ignite
            (seconds).
    """

    def __init__(self, spot_delay_s: float, limits: tuple):
        """Initialize the Perryman spotting model.

        Args:
            spot_delay_s (float): Spot fire ignition delay (seconds).
            limits (tuple): ``(x_lim, y_lim)`` simulation domain bounds
                (meters).
        """
        self.embers = []
        self.spot_delay_s = spot_delay_s

        self.x_lim, self.y_lim = limits

    def loft(self, cell: Cell):
        """Sample firebrand landing locations from a burning cell.

        Compute downwind (parallel) and crosswind (perpendicular) landing
        distances, combine them into absolute coordinates, and store
        firebrands that land beyond 50 meters from the source.

        Args:
            cell (Cell): Burning cell to generate firebrands from.

        Side Effects:
            Sets ``cell.lofted = True``. Appends ember dicts to
            ``self.embers``.
        """
        cell.lofted = True

        wind_speed, wind_dir_deg = cell.curr_wind()

        if wind_speed == 0:
            return

        # Compute distances parallel and perpendicular to wind dir
        par_distances = self.compute_parallel_distances(cell)
        perp_distances = self.compute_perp_distances(cell)

        # Generate actual coordinate where each brand will land

        # Get the downwind distances
        wind_dir = np.deg2rad(wind_dir_deg)
        down_wind_x_vec = par_distances * np.sin(wind_dir)
        down_wind_y_vec = par_distances * np.cos(wind_dir)

        # Get the distances perpendicular to wind direction
        perp_dir = wind_dir + (np.pi / 2)
        perp_x_vec = perp_distances * np.sin(perp_dir)
        perp_y_vec = perp_distances * np.cos(perp_dir)

        # Add the vectors
        tot_vec_x = down_wind_x_vec + perp_x_vec
        tot_vec_y = down_wind_y_vec + perp_y_vec

        # Compute absolute coordinates relative to the cell position
        x_coords = cell.x_pos + tot_vec_x
        y_coords = cell.y_pos + tot_vec_y

        # Ensure coordinates are within the limits
        x_coords = np.clip(x_coords, 0, self.x_lim)
        y_coords = np.clip(y_coords, 0, self.y_lim)

        # Calculate the distance from the cell center to each ember
        d = np.sqrt((cell.x_pos - x_coords)**2 + (cell.y_pos - y_coords)**2)

        # Filter out distances below minimum spotting distance
        filtered_indices = np.where(d > 50)[0]
        filtered_d = d[filtered_indices]
        x_coords = x_coords[filtered_indices]
        y_coords = y_coords[filtered_indices]

        for i in range(len(filtered_d)):
            ember = {
                'x': x_coords[i],
                'y': y_coords[i],
                'd': filtered_d[i]
            }

            self.embers.append(ember)

    def compute_parallel_distances(self, cell: Cell) -> np.ndarray:
        """Sample downwind (parallel) firebrand distances.

        Compute the canopy-height wind speed, fireline intensity, and
        Froude number, then sample distances from a lognormal distribution
        whose parameters depend on the Froude number regime.

        Args:
            cell (Cell): Burning cell providing wind, canopy, and intensity.

        Returns:
            np.ndarray: Downwind distances for each firebrand (meters),
                shape ``(num_firebrands,)``.
        """
        wind_spd_ft_s = m_to_ft(cell.curr_wind()[0])
        canopy_ht = m_to_ft(cell.canopy_height)
        u_h_ft_s = wind_spd_ft_s / (np.log((20 + 0.36 * canopy_ht)/(0.1313 * canopy_ht)))
        u_h = ft_to_m(u_h_ft_s)

        # Get the fireline intensity in kW/m
        I_f_btu_ft_min = np.max(cell.I_ss)
        I_f_kW_m = BTU_ft_min_to_kW_m(I_f_btu_ft_min)

        # Compute Froude number to compute the lognormal distribution
        denom = np.sqrt(g * (I_f_kW_m / (rN * cpg * TN * np.sqrt(g)))**(2/3.0))

        if denom == 0:
            Fr = float('inf')
        else:
            Fr = u_h / denom

        # Set lognormal distribution parameters
        if Fr <= 1:
            sigma = 0.86 * (I_f_kW_m**(-0.21)*u_h**(0.44)) + 0.19
            mu =    1.47 * (I_f_kW_m**(0.54)*u_h**(-0.55)) + 1.14

        else:
            sigma = 4.95 * (I_f_kW_m**(-0.01)*u_h**(-0.02)) - 3.48
            mu =    1.32 * (I_f_kW_m**(0.26)*u_h**(0.11)) - 0.02

        # Sample distances for firebrands
        parallel_distances = np.random.lognormal(mean=mu, sigma=sigma, size=num_firebrands)

        return parallel_distances

    def compute_perp_distances(self, cell: Cell) -> np.ndarray:
        """Sample crosswind (perpendicular) firebrand distances.

        Use a normal distribution centered at zero with standard deviation
        equal to half the cell's flat-to-flat width.

        Args:
            cell (Cell): Cell providing ``cell_size`` for scale computation.

        Returns:
            np.ndarray: Crosswind distances for each firebrand (meters),
                shape ``(num_firebrands,)``.
        """
        # Horizontal standard dev. is half of cell width
        flat_to_flat = cell.cell_size * (np.sqrt(3)/2)
        sv = flat_to_flat / 2

        perp_distances = np.random.normal(loc=0, scale=sv, size=num_firebrands)

        return perp_distances

__init__(spot_delay_s, limits)

Initialize the Perryman spotting model.

Parameters:

Name Type Description Default
spot_delay_s float

Spot fire ignition delay (seconds).

required
limits tuple

(x_lim, y_lim) simulation domain bounds (meters).

required
Source code in embrs/models/perryman_spot.py
42
43
44
45
46
47
48
49
50
51
52
53
def __init__(self, spot_delay_s: float, limits: tuple):
    """Initialize the Perryman spotting model.

    Args:
        spot_delay_s (float): Spot fire ignition delay (seconds).
        limits (tuple): ``(x_lim, y_lim)`` simulation domain bounds
            (meters).
    """
    self.embers = []
    self.spot_delay_s = spot_delay_s

    self.x_lim, self.y_lim = limits

compute_parallel_distances(cell)

Sample downwind (parallel) firebrand distances.

Compute the canopy-height wind speed, fireline intensity, and Froude number, then sample distances from a lognormal distribution whose parameters depend on the Froude number regime.

Parameters:

Name Type Description Default
cell Cell

Burning cell providing wind, canopy, and intensity.

required

Returns:

Type Description
ndarray

np.ndarray: Downwind distances for each firebrand (meters), shape (num_firebrands,).

Source code in embrs/models/perryman_spot.py
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
def compute_parallel_distances(self, cell: Cell) -> np.ndarray:
    """Sample downwind (parallel) firebrand distances.

    Compute the canopy-height wind speed, fireline intensity, and
    Froude number, then sample distances from a lognormal distribution
    whose parameters depend on the Froude number regime.

    Args:
        cell (Cell): Burning cell providing wind, canopy, and intensity.

    Returns:
        np.ndarray: Downwind distances for each firebrand (meters),
            shape ``(num_firebrands,)``.
    """
    wind_spd_ft_s = m_to_ft(cell.curr_wind()[0])
    canopy_ht = m_to_ft(cell.canopy_height)
    u_h_ft_s = wind_spd_ft_s / (np.log((20 + 0.36 * canopy_ht)/(0.1313 * canopy_ht)))
    u_h = ft_to_m(u_h_ft_s)

    # Get the fireline intensity in kW/m
    I_f_btu_ft_min = np.max(cell.I_ss)
    I_f_kW_m = BTU_ft_min_to_kW_m(I_f_btu_ft_min)

    # Compute Froude number to compute the lognormal distribution
    denom = np.sqrt(g * (I_f_kW_m / (rN * cpg * TN * np.sqrt(g)))**(2/3.0))

    if denom == 0:
        Fr = float('inf')
    else:
        Fr = u_h / denom

    # Set lognormal distribution parameters
    if Fr <= 1:
        sigma = 0.86 * (I_f_kW_m**(-0.21)*u_h**(0.44)) + 0.19
        mu =    1.47 * (I_f_kW_m**(0.54)*u_h**(-0.55)) + 1.14

    else:
        sigma = 4.95 * (I_f_kW_m**(-0.01)*u_h**(-0.02)) - 3.48
        mu =    1.32 * (I_f_kW_m**(0.26)*u_h**(0.11)) - 0.02

    # Sample distances for firebrands
    parallel_distances = np.random.lognormal(mean=mu, sigma=sigma, size=num_firebrands)

    return parallel_distances

compute_perp_distances(cell)

Sample crosswind (perpendicular) firebrand distances.

Use a normal distribution centered at zero with standard deviation equal to half the cell's flat-to-flat width.

Parameters:

Name Type Description Default
cell Cell

Cell providing cell_size for scale computation.

required

Returns:

Type Description
ndarray

np.ndarray: Crosswind distances for each firebrand (meters), shape (num_firebrands,).

Source code in embrs/models/perryman_spot.py
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
def compute_perp_distances(self, cell: Cell) -> np.ndarray:
    """Sample crosswind (perpendicular) firebrand distances.

    Use a normal distribution centered at zero with standard deviation
    equal to half the cell's flat-to-flat width.

    Args:
        cell (Cell): Cell providing ``cell_size`` for scale computation.

    Returns:
        np.ndarray: Crosswind distances for each firebrand (meters),
            shape ``(num_firebrands,)``.
    """
    # Horizontal standard dev. is half of cell width
    flat_to_flat = cell.cell_size * (np.sqrt(3)/2)
    sv = flat_to_flat / 2

    perp_distances = np.random.normal(loc=0, scale=sv, size=num_firebrands)

    return perp_distances

loft(cell)

Sample firebrand landing locations from a burning cell.

Compute downwind (parallel) and crosswind (perpendicular) landing distances, combine them into absolute coordinates, and store firebrands that land beyond 50 meters from the source.

Parameters:

Name Type Description Default
cell Cell

Burning cell to generate firebrands from.

required
Side Effects

Sets cell.lofted = True. Appends ember dicts to self.embers.

Source code in embrs/models/perryman_spot.py
 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
def loft(self, cell: Cell):
    """Sample firebrand landing locations from a burning cell.

    Compute downwind (parallel) and crosswind (perpendicular) landing
    distances, combine them into absolute coordinates, and store
    firebrands that land beyond 50 meters from the source.

    Args:
        cell (Cell): Burning cell to generate firebrands from.

    Side Effects:
        Sets ``cell.lofted = True``. Appends ember dicts to
        ``self.embers``.
    """
    cell.lofted = True

    wind_speed, wind_dir_deg = cell.curr_wind()

    if wind_speed == 0:
        return

    # Compute distances parallel and perpendicular to wind dir
    par_distances = self.compute_parallel_distances(cell)
    perp_distances = self.compute_perp_distances(cell)

    # Generate actual coordinate where each brand will land

    # Get the downwind distances
    wind_dir = np.deg2rad(wind_dir_deg)
    down_wind_x_vec = par_distances * np.sin(wind_dir)
    down_wind_y_vec = par_distances * np.cos(wind_dir)

    # Get the distances perpendicular to wind direction
    perp_dir = wind_dir + (np.pi / 2)
    perp_x_vec = perp_distances * np.sin(perp_dir)
    perp_y_vec = perp_distances * np.cos(perp_dir)

    # Add the vectors
    tot_vec_x = down_wind_x_vec + perp_x_vec
    tot_vec_y = down_wind_y_vec + perp_y_vec

    # Compute absolute coordinates relative to the cell position
    x_coords = cell.x_pos + tot_vec_x
    y_coords = cell.y_pos + tot_vec_y

    # Ensure coordinates are within the limits
    x_coords = np.clip(x_coords, 0, self.x_lim)
    y_coords = np.clip(y_coords, 0, self.y_lim)

    # Calculate the distance from the cell center to each ember
    d = np.sqrt((cell.x_pos - x_coords)**2 + (cell.y_pos - y_coords)**2)

    # Filter out distances below minimum spotting distance
    filtered_indices = np.where(d > 50)[0]
    filtered_d = d[filtered_indices]
    x_coords = x_coords[filtered_indices]
    y_coords = y_coords[filtered_indices]

    for i in range(len(filtered_d)):
        ember = {
            'x': x_coords[i],
            'y': y_coords[i],
            'd': filtered_d[i]
        }

        self.embers.append(ember)

BTU_ft2_min_to_kW_m2(f_btu_ft2_min)

Convert heat flux from BTU/(ft^2*min) to kW/m^2.

Parameters:

Name Type Description Default
f_btu_ft2_min float

Heat flux in BTU/(ft^2*min).

required

Returns:

Name Type Description
float float

Heat flux in kW/m^2.

Source code in embrs/utilities/unit_conversions.py
220
221
222
223
224
225
226
227
228
229
def BTU_ft2_min_to_kW_m2(f_btu_ft2_min: float) -> float:
    """Convert heat flux from BTU/(ft^2*min) to kW/m^2.

    Args:
        f_btu_ft2_min (float): Heat flux in BTU/(ft^2*min).

    Returns:
        float: Heat flux in kW/m^2.
    """
    return f_btu_ft2_min * _BTU_FT2_MIN_TO_KW_M2

BTU_ft_min_to_kW_m(f_btu_ft_min)

Convert fireline intensity from BTU/(ft*min) to kW/m.

Parameters:

Name Type Description Default
f_btu_ft_min float

Fireline intensity in BTU/(ft*min).

required

Returns:

Name Type Description
float float

Fireline intensity in kW/m.

Source code in embrs/utilities/unit_conversions.py
232
233
234
235
236
237
238
239
240
241
def BTU_ft_min_to_kW_m(f_btu_ft_min: float) -> float:
    """Convert fireline intensity from BTU/(ft*min) to kW/m.

    Args:
        f_btu_ft_min (float): Fireline intensity in BTU/(ft*min).

    Returns:
        float: Fireline intensity in kW/m.
    """
    return f_btu_ft_min * _BTU_FT_MIN_TO_KW_M

BTU_ft_min_to_kcal_s_m(f_btu_ft_min)

Convert fireline intensity from BTU/(ftmin) to kcal/(sm).

Parameters:

Name Type Description Default
f_btu_ft_min float

Fireline intensity in BTU/(ft*min).

required

Returns:

Name Type Description
float float

Fireline intensity in kcal/(s*m).

Source code in embrs/utilities/unit_conversions.py
244
245
246
247
248
249
250
251
252
253
def BTU_ft_min_to_kcal_s_m(f_btu_ft_min: float) -> float:
    """Convert fireline intensity from BTU/(ft*min) to kcal/(s*m).

    Args:
        f_btu_ft_min (float): Fireline intensity in BTU/(ft*min).

    Returns:
        float: Fireline intensity in kcal/(s*m).
    """
    return f_btu_ft_min * _BTU_FT_MIN_TO_KCAL_S_M

BTU_lb_to_cal_g(f_btu_lb)

Convert heat content from BTU/lb to cal/g.

Parameters:

Name Type Description Default
f_btu_lb float

Heat content in BTU/lb.

required

Returns:

Name Type Description
float float

Heat content in cal/g.

Source code in embrs/utilities/unit_conversions.py
272
273
274
275
276
277
278
279
280
281
def BTU_lb_to_cal_g(f_btu_lb: float) -> float:
    """Convert heat content from BTU/lb to cal/g.

    Args:
        f_btu_lb (float): Heat content in BTU/lb.

    Returns:
        float: Heat content in cal/g.
    """
    return f_btu_lb * _BTU_LB_TO_CAL_G

F_to_C(f_f)

Convert temperature from Fahrenheit to Celsius.

Parameters:

Name Type Description Default
f_f float

Temperature in degrees Fahrenheit.

required

Returns:

Name Type Description
float float

Temperature in degrees Celsius.

Source code in embrs/utilities/unit_conversions.py
48
49
50
51
52
53
54
55
56
57
def F_to_C(f_f: float) -> float:
    """Convert temperature from Fahrenheit to Celsius.

    Args:
        f_f (float): Temperature in degrees Fahrenheit.

    Returns:
        float: Temperature in degrees Celsius.
    """
    return (5.0 / 9.0) * (f_f - 32.0)

KiSq_to_Lbsft2(f_kisq)

Convert fuel loading from kg/m^2 to lb/ft^2.

Parameters:

Name Type Description Default
f_kisq float

Fuel loading in kg/m^2.

required

Returns:

Name Type Description
float float

Fuel loading in lb/ft^2.

Source code in embrs/utilities/unit_conversions.py
156
157
158
159
160
161
162
163
164
165
def KiSq_to_Lbsft2(f_kisq: float) -> float:
    """Convert fuel loading from kg/m^2 to lb/ft^2.

    Args:
        f_kisq (float): Fuel loading in kg/m^2.

    Returns:
        float: Fuel loading in lb/ft^2.
    """
    return f_kisq * _KISQ_TO_LBSFT2

KiSq_to_TPA(f_kisq)

Convert fuel loading from kg/m^2 to tons per acre.

Parameters:

Name Type Description Default
f_kisq float

Fuel loading in kg/m^2.

required

Returns:

Name Type Description
float float

Fuel loading in tons per acre.

Source code in embrs/utilities/unit_conversions.py
204
205
206
207
208
209
210
211
212
213
def KiSq_to_TPA(f_kisq: float) -> float:
    """Convert fuel loading from kg/m^2 to tons per acre.

    Args:
        f_kisq (float): Fuel loading in kg/m^2.

    Returns:
        float: Fuel loading in tons per acre.
    """
    return f_kisq * _KISQ_TO_TPA

Lbsft2_to_KiSq(f_libsft2)

Convert fuel loading from lb/ft^2 to kg/m^2.

Parameters:

Name Type Description Default
f_libsft2 float

Fuel loading in lb/ft^2.

required

Returns:

Name Type Description
float float

Fuel loading in kg/m^2.

Source code in embrs/utilities/unit_conversions.py
144
145
146
147
148
149
150
151
152
153
def Lbsft2_to_KiSq(f_libsft2: float) -> float:
    """Convert fuel loading from lb/ft^2 to kg/m^2.

    Args:
        f_libsft2 (float): Fuel loading in lb/ft^2.

    Returns:
        float: Fuel loading in kg/m^2.
    """
    return f_libsft2 * _LBSFT2_TO_KISQ

Lbsft2_to_TPA(f_lbsft2)

Convert fuel loading from lb/ft^2 to tons per acre.

Parameters:

Name Type Description Default
f_lbsft2 float

Fuel loading in lb/ft^2.

required

Returns:

Name Type Description
float float

Fuel loading in tons per acre.

Source code in embrs/utilities/unit_conversions.py
192
193
194
195
196
197
198
199
200
201
def Lbsft2_to_TPA(f_lbsft2: float) -> float:
    """Convert fuel loading from lb/ft^2 to tons per acre.

    Args:
        f_lbsft2 (float): Fuel loading in lb/ft^2.

    Returns:
        float: Fuel loading in tons per acre.
    """
    return f_lbsft2 * _LBSFT2_TO_TPA

TPA_to_KiSq(f_tpa)

Convert fuel loading from tons per acre to kg/m^2.

Parameters:

Name Type Description Default
f_tpa float

Fuel loading in tons per acre.

required

Returns:

Name Type Description
float float

Fuel loading in kg/m^2.

Source code in embrs/utilities/unit_conversions.py
168
169
170
171
172
173
174
175
176
177
def TPA_to_KiSq(f_tpa: float) -> float:
    """Convert fuel loading from tons per acre to kg/m^2.

    Args:
        f_tpa (float): Fuel loading in tons per acre.

    Returns:
        float: Fuel loading in kg/m^2.
    """
    return f_tpa / _TPA_TO_KISQ_DIVISOR

TPA_to_Lbsft2(f_tpa)

Convert fuel loading from tons per acre to lb/ft^2.

Parameters:

Name Type Description Default
f_tpa float

Fuel loading in tons per acre.

required

Returns:

Name Type Description
float float

Fuel loading in lb/ft^2.

Source code in embrs/utilities/unit_conversions.py
180
181
182
183
184
185
186
187
188
189
def TPA_to_Lbsft2(f_tpa: float) -> float:
    """Convert fuel loading from tons per acre to lb/ft^2.

    Args:
        f_tpa (float): Fuel loading in tons per acre.

    Returns:
        float: Fuel loading in lb/ft^2.
    """
    return f_tpa * _TPA_TO_LBSFT2

cal_g_to_BTU_lb(f_cal_g)

Convert heat content from cal/g to BTU/lb.

Parameters:

Name Type Description Default
f_cal_g float

Heat content in cal/g.

required

Returns:

Name Type Description
float float

Heat content in BTU/lb.

Source code in embrs/utilities/unit_conversions.py
260
261
262
263
264
265
266
267
268
269
def cal_g_to_BTU_lb(f_cal_g: float) -> float:
    """Convert heat content from cal/g to BTU/lb.

    Args:
        f_cal_g (float): Heat content in cal/g.

    Returns:
        float: Heat content in BTU/lb.
    """
    return f_cal_g * _CAL_G_TO_BTU_LB

ft_min_to_m_s(f_ft_min)

Convert speed from feet per minute to meters per second.

Parameters:

Name Type Description Default
f_ft_min float

Speed in ft/min.

required

Returns:

Name Type Description
float float

Speed in m/s.

Source code in embrs/utilities/unit_conversions.py
 92
 93
 94
 95
 96
 97
 98
 99
100
101
def ft_min_to_m_s(f_ft_min: float) -> float:
    """Convert speed from feet per minute to meters per second.

    Args:
        f_ft_min (float): Speed in ft/min.

    Returns:
        float: Speed in m/s.
    """
    return f_ft_min * _FT_MIN_TO_M_S

ft_min_to_mph(f_ft_min)

Convert speed from feet per minute to miles per hour.

Parameters:

Name Type Description Default
f_ft_min float

Speed in ft/min.

required

Returns:

Name Type Description
float float

Speed in mph.

Source code in embrs/utilities/unit_conversions.py
116
117
118
119
120
121
122
123
124
125
def ft_min_to_mph(f_ft_min: float) -> float:
    """Convert speed from feet per minute to miles per hour.

    Args:
        f_ft_min (float): Speed in ft/min.

    Returns:
        float: Speed in mph.
    """
    return f_ft_min * _FT_MIN_TO_MPH

ft_to_m(f_ft)

Convert length from feet to meters.

Parameters:

Name Type Description Default
f_ft float

Length in feet.

required

Returns:

Name Type Description
float float

Length in meters.

Source code in embrs/utilities/unit_conversions.py
76
77
78
79
80
81
82
83
84
85
def ft_to_m(f_ft: float) -> float:
    """Convert length from feet to meters.

    Args:
        f_ft (float): Length in feet.

    Returns:
        float: Length in meters.
    """
    return f_ft * _FT_TO_M

m_s_to_ft_min(m_s)

Convert speed from meters per second to feet per minute.

Parameters:

Name Type Description Default
m_s float

Speed in m/s.

required

Returns:

Name Type Description
float float

Speed in ft/min.

Source code in embrs/utilities/unit_conversions.py
104
105
106
107
108
109
110
111
112
113
def m_s_to_ft_min(m_s: float) -> float:
    """Convert speed from meters per second to feet per minute.

    Args:
        m_s (float): Speed in m/s.

    Returns:
        float: Speed in ft/min.
    """
    return m_s * _M_S_TO_FT_MIN

m_to_ft(f_m)

Convert length from meters to feet.

Parameters:

Name Type Description Default
f_m float

Length in meters.

required

Returns:

Name Type Description
float float

Length in feet.

Source code in embrs/utilities/unit_conversions.py
64
65
66
67
68
69
70
71
72
73
def m_to_ft(f_m: float) -> float:
    """Convert length from meters to feet.

    Args:
        f_m (float): Length in meters.

    Returns:
        float: Length in feet.
    """
    return f_m * _M_TO_FT

mph_to_ft_min(f_mph)

Convert speed from miles per hour to feet per minute.

Parameters:

Name Type Description Default
f_mph float

Speed in mph.

required

Returns:

Name Type Description
float float

Speed in ft/min.

Source code in embrs/utilities/unit_conversions.py
128
129
130
131
132
133
134
135
136
137
def mph_to_ft_min(f_mph: float) -> float:
    """Convert speed from miles per hour to feet per minute.

    Args:
        f_mph (float): Speed in mph.

    Returns:
        float: Speed in ft/min.
    """
    return f_mph * _MPH_TO_FT_MIN