-
Notifications
You must be signed in to change notification settings - Fork 18
/
module_fr_sfire_core.F
1232 lines (1091 loc) · 40.7 KB
/
module_fr_sfire_core.F
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!
!*** Jan Mandel August-October 2007 email: jmandel@ucar.edu or Jan.Mandel@gmail.com
!
! With contributions by Minjeong Kim.
!#define DEBUG_OUT
module module_fr_sfire_core
use module_fr_sfire_phys
use module_fr_sfire_util
! The mathematical core of the fire spread model. No physical constants here.
!
! subroutine sfire_core: only this routine should be called from the outside.
! subroutine fuel_left: compute remaining fuel from time of ignition.
! subroutine prop_ls: propagation of curve in normal direction.
contains
!
!****************************************
!
subroutine init_no_fire(&
ifds,ifde,jfds,jfde, &
ifms,ifme,jfms,jfme, &
ifts,ifte,jfts,jfte, &
fdx,fdy,time_now, & ! scalars in
fuel_frac,lfn,tign) ! arrays out
implicit none
!*** purpose: initialize model to no fire
!*** arguments
integer, intent(in):: ifds,ifde,jfds,jfde ! fire domain bounds
integer, intent(in):: ifts,ifte,jfts,jfte ! fire tile bounds
integer, intent(in):: ifms,ifme,jfms,jfme ! array bounds
real, intent(in) :: fdx,fdy,time_now ! mesh spacing, time
real, intent(out), dimension (ifms:ifme,jfms:jfme) :: &
fuel_frac,lfn,tign ! model state
!*** calls
intrinsic epsilon
!*** local
integer:: i,j
real lfn_init,time_init
do j=jfts,jfte
do i=ifts,ifte
fuel_frac(i,j)=1. ! fuel at start is 1 by definition
enddo
enddo
lfn_init = 2*max((ifde-ifds+1)*fdx,(jfde-jfds+1)*fdy) ! more than domain diameter
time_init=time_now + max(time_now,1.0)*epsilon(time_now) ! a bit in future
do j=jfts,jfte
do i=ifts,ifte
tign(i,j) = time_init ! ignition in future
lfn(i,j) = lfn_init ! no fire
enddo
enddo
call message('init_model_no_fire: state set to no fire')
end subroutine init_no_fire
!
!******************
!
subroutine ignite_fire( ifds,ifde,jfds,jfde, & ! fire domain dims - the whole domain
ifms,ifme,jfms,jfme, &
ifts,ifte,jfts,jfte, &
sx,sy,ex,ey,r,time_ign,fdx,fdy, &
lfn,tign,ignited)
implicit none
!*** purpose: ignite a circular/line fire
!*** arguments
integer, intent(in):: ifds,ifde,jfds,jfde ! fire domain bounds
integer, intent(in):: ifts,ifte,jfts,jfte ! fire tile bounds
integer, intent(in):: ifms,ifme,jfms,jfme ! array bounds
real, intent(in):: time_ign ! the ignition time of the fire
real, intent(in):: sx,sy ! start of ignition line, from lower left corner
real, intent(in):: ex,ey ! end of ignition line, or zero
real, intent(in):: r ! all within the radius of the line will ignite
real, intent(in):: fdx,fdy ! mesh spacing (m)
real, intent(inout), dimension (ifms:ifme,jfms:jfme) :: &
lfn, tign ! level function, ignition time (state)
integer, intent(out):: ignited ! number of nodes newly ignited
!*** local
integer:: i,j
real::mx,my,ax,ay,dam2,d,dames,des2,am_es,cos2,lfn_new,dmc2
logical::point
character(len=128):: msg
ignited=0
point = ex .eq. 0.0 .or. ey .eq. 0.0
if (.not.point)then
! midpoint m = (mx,my)
mx = (sx + ex)/2
my = (sy + ey)/2
else
mx = sx
my = sy
endif
do j=jfts,jfte
do i=ifts,ifte
! coordinates of the node (i,j), the lower left corner of the domain is (0 0)
ax = fdx*(i - ifds + 0.5)
ay = fdy*(j - jfds + 0.5)
dam2=(ax-mx)*(ax-mx)+(ay-my)*(ay-my) ! |a-m|^2
if(point)then
d=sqrt(dam2)
else
! compute distance as distance from midpoint minus correction
! |a-c|^2 = |a-m|^2 - |m-c|^2
! when |m-c| >= |s-e|/2 use distance from the endpoint instead
!
! a
! /| \
! s---m-c--e
!
! |m-c| = |a-m| cos (a-m,e-s)
! = |a-m| (a-m).(e-s))/(|a-m|*|e-s|)
des2 = (ex-sx)*(ex-sx)+(ey-sy)*(ey-sy) ! |e-s|^2
dames = dam2*des2
if(dames>0)then
am_es=(ax-mx)*(ex-sx)+(ay-my)*(ey-sy) ! (a-m).(e-s)
cos2 = (am_es*am_es)/dames ! cos^2 (a-m,e-s)
else
cos2 = 0.
endif
dmc2 = dam2*cos2 ! |m-c|^2
if(4.*dmc2 <= des2)then
d = sqrt(max(dam2 - dmc2,0.)) ! just in case, rounding
elseif(am_es>0)then ! cos > 0, closest is e
d = sqrt((ax-ex)*(ax-ex)+(ay-ey)*(ay-ey)) ! |a-e|
else
d = sqrt((ax-sx)*(ax-sx)+(ay-sy)*(ay-sy)) ! |a-s|
endif
endif
lfn_new=d-r
if(lfn(i,j)>0 .and. lfn_new<=0) then
tign(i,j)=time_ign ! newly ignited now
ignited=ignited+1 ! count
endif
lfn(i,j)=min(lfn(i,j),lfn_new) ! update the level set function
enddo
enddo
write(msg,'(a,2f10.1,a,2f10.1,a,f8.1,a,f8.1,a,i6)')'ignite_fire: from',sx,sy,' to ',&
ex,ey,' radius ',r,' time',time_ign,' ignited nodes',ignited
call message(msg)
end subroutine ignite_fire
!
!**********************
!
subroutine fuel_left( &
ims,ime,jms,jme, &
its,ite,jts,jte, &
ifs,ife,jfs,jfe, &
ir,jr, &
lfn, tign, fuel_time, time_now, fuel_frac)
implicit none
!*** purpose: determine fraction of fuel remaining
!*** Jan Mandel August 2007 email: jmandel@ucar.edu or Jan.Mandel@gmail.com
!*** arguments
integer, intent(in) :: its,ite,jts,jte,ims,ime,jms,jme,ifs,ife,jfs,jfe,ir,jr
real, intent(in), dimension(ims:ime,jms:jme)::lfn,tign,fuel_time
real, intent(in):: time_now
real, intent(out), dimension(ifs:ife,jfs:jfe)::fuel_frac
! ims,ime,jms,jme in memory dimensions
! its,ite,jts,jte in tile dimensions (cells where fuel_frac computed)
! ifs,ife,jfs,jfe in fuel_frac memory dimensions
! ir,jr, in refinement - quadrature cells per fire cell, even, at least 2
! lfn in level function, at nodes at midpoints of cells
! tign in ignition time, at nodes at nodes at midpoints of cells
! fuel_time in time constant of fuel, per cell
! time_now in time now
! fuel_frac out fraction of fuel remaining, per cell
!*** local
integer::i,j,iccell,jccell
real::xfcell,yfcell
#define IFCELLS (ite-its+1)*ir
#define JFCELLS (jte-jts+1)*jr
real,dimension(0:IFCELLS+1,0:JFCELLS+1)::lff,tif
real,dimension(1:IFCELLS,1:JFCELLS)::fuel_left_f
call check_mesh_2dim(its-1,ite+1,jts-1,jte+1,ims,ime,jms,jme)
call check_mesh_2dim(its,ite,jts,jte,ifs,ife,jfs,jfe)
! interpolate level set function to finer grid
call interpolate_2d( &
ims,ime,jms,jme, & ! array coarse grid
its-1,ite+1,jts-1,jte+1, & ! dimensions coarse grid
0,IFCELLS+1,0,JFCELLS+1, & ! array fine grid
0,IFCELLS+1,0,JFCELLS+1, & ! dimensions fine grid
ir,jr, & ! refinement ratio
real(its),real(jts),ir/2.,jr/2., & ! (ip2,jp2) on coarse lines up with (ip1,jp1) on fine grid
lfn, & ! in coarse grid
lff ) ! out fine grid
call interpolate_2d( &
ims,ime,jms,jme, & ! array coarse grid
its-1,ite+1,jts-1,jte+1, & ! dimensions coarse grid
0,IFCELLS+1,0,JFCELLS+1, & ! array fine grid
0,IFCELLS+1,0,JFCELLS+1, & ! dimensions fine grid
ir,jr, & ! refinement ratio
real(its),real(jts),ir/2.,jr/2., & ! (ip2,jp2) on coarse lines up with (ip1,jp1) on fine grid
tign, & ! in coarse grid
tif ) ! out fine grid
!
! example for ir=2:
!
! its-1 coarse cell its ite ite+1
! -------X------------|-------------X-------------|--........----|--------X-------|--------X
! fine node 0 1 2 2*(ite-its)
! fine cell 1 2 2*(ite-its+1)
do j=1,JFCELLS
do i=1,IFCELLS
xfcell=i-0.5 ! middle of the fine cell in fine grid coordinates
yfcell=j-0.5
iccell=nint(its+xfcell/ir-0.5) ! the indices of the nearest coarse midpoint
jccell=nint(jts+yfcell/jr-0.5) ! use fuel_time of the coarse cell the fine cell is in
fuel_left_f(i,j)=fuel_left_cell( &
lff(i,j),lff(i,j+1),lff(i+1,j),lff(i+1,j+1),&
tif(i,j),tif(i,j+1),tif(i+1,j),tif(i+1,j+1),&
time_now, fuel_time(iccell,jccell)) &
/(ir*jr) ! premultiply by weight
enddo
enddo
call sum_2d_cells( &
1,IFCELLS,1,JFCELLS, &
1,IFCELLS,1,JFCELLS, &
fuel_left_f, & ! input
ifs,ife,jfs,jfe, &
its,ite,jts,jte, &
fuel_frac) ! output
end subroutine fuel_left
!
!************************
!
real function fuel_left_cell( &
lfn00,lfn01,lfn10,lfn11, &
tign00,tign01,tign10,tign11,&
time_now, fuel_time_cell)
!*** purpose: compute the fuel fraction left in one cell
implicit none
!*** arguments
real, intent(in)::lfn00,lfn01,lfn10,lfn11 ! level set function at 4 corners of the cell
real, intent(in)::tign00,tign01,tign10,tign11! ignition time at the 4 corners of the cell
real, intent(in)::time_now ! the time now
real, intent(in)::fuel_time_cell ! time to burns off to 1/e
!*** Description
! The area burning is given by the condition L <= 0, where the function P is
! interpolated from lfn(i,j)
!
! The time since ignition is the function T, interpolated in from tign(i,j)-time_now.
! The values of tign(i,j) where lfn(i,j)>=0 are ignored, tign(i,j)=0 is taken
! when lfn(i,j)=0.
!
! The function computes an approxmation of the integral
!
!
! /\
! |
! fuel_frac_left = 1 - | 1 - exp(-T(x,y)*fuel_speed)) dxdy
! |
! \/
! 0<x<1
! 0<y<1
! L(x,y)<=0
!
! When the cell is not burning at all (all lfn>=0), then fuel_frac(i,j)=1.
! Because of symmetries, the result should not depend on the mesh spacing dx dy
! so dx=1 and dy=1 assumed.
!
! Example:
!
! lfn<0 lfn>0
! (0,1) -----O--(1,1) O = points on the fireline, T=tnow
! | \ | A = the burning area for computing
! | \| fuel_frac(i,j)
! | A O
! | |
! | |
! (0,0)---------(1,0)
! lfn<0 lfn<0
!
! Approximations allowed:
! The fireline can be approximated by straight line(s).
! When all cell is burning, approximation by 1 point Gaussian quadrature is OK.
!
! Requirements:
! 1. The output should be a continuous function of the arrays lfn and
! tign whenever lfn(i,j)=0 implies tign(i,j)=tnow.
! 2. The output should be invariant to the symmetries of the input in each cell.
! 3. Arbitrary combinations of the signs of lfn(i,j) should work.
! 4. The result should be at least 1st order accurate in the sense that it is
! exact if the time from ignition is a linear function.
!
! If time from ignition is approximated by polynomial in the burnt
! region of the cell, this is integral of polynomial times exponential
! over a polygon, which can be computed exactly.
!
! Requirement 4 is particularly important when there is a significant decrease
! of the fuel fraction behind the fireline on the mesh scale, because the
! rate of fuel decrease right behind the fireline is much larger
! (exponential...). This will happen when
!
! change of time from ignition within one mesh cell * fuel speed is not << 1
!
! This is the same as
!
! mesh cell size*fuel_speed
! ------------------------- is not << 1
! fireline speed
!
!
! When X is large then the fuel burnt in one timestep in the cell is
! approximately proportional to length of fireline in that cell.
!
! When X is small then the fuel burnt in one timestep in the cell is
! approximately proportional to the area of the burning region.
!
!*** calls
intrinsic tiny
!*** local
real::ps,aps,area,ta,out
real::t00,t01,t10,t11
real,parameter::safe=tiny(aps)
character(len=128)::msg
!*** local
integer::i,j,k
! least squares
integer::mmax,nb,nmax,pmax,nin,nout
parameter(mmax=3,nb=64,nmax=8,pmax=8)
integer lda, ldb, lwork, info
parameter (lda=nmax, ldb=nmax, lwork=mmax+nmax+nb*(nmax+pmax))
integer n,m,p
real,dimension(lda,mmax):: mA
real,dimension(nmax):: vecD
real,dimension(lwork):: WORK
real,dimension(pmax):: vecY
real,dimension(mmax):: vecX
real,dimension(ldb,pmax)::mB
real,dimension(mmax)::u
real::tweight,tdist
integer::kk,ll,ss
real::rnorm
real,dimension(8,2)::xylist,xytlist
real,dimension(8)::tlist,llist,xt
real,dimension(5)::xx,yy
real,dimension(5)::lfn,tign
integer:: npoint
real::tt,x0,y0,xts,xte,yts,yte,xt1,xt2
real::lfn0,lfn1,dist,nr,c,s,errQ,ae,ce,ceae,a0,a1,a2,d,cet
real::s1,s2,s3
real::upper,lower,ah,ch,aa,cc,aupp,cupp,alow,clow
real,dimension(2,2)::mQ
real,dimension(2)::ut
!calls
intrinsic epsilon
real, parameter:: zero=0.,one=1.,eps=epsilon(zero)
! external functions
real::snrm2
double precision :: dnrm2
external dnrm2
external snrm2
! external subroutines
external sggglm
external dggglm
!executable statements
! a very crude approximation - replace by a better code
! call check_mesh_2dim(ids,ide+1,jds,jde+1,ims,ime,jms,jme)
t00=tign00-time_now
if(lfn00>=0. .or. t00>0.)t00=0.
t01=tign01-time_now
if(lfn01>=0. .or. t01>0.)t01=0.
t10=tign10-time_now
if(lfn10>=0. .or. t10>0.)t10=0.
t11=tign11-time_now
if(lfn11>=0. .or. t11>0.)t11=0.
!if (t00<0 .or. t01 <0 .or. t10<0 .or. t11<0) then
! print *,'tign00=',tign00,'tign10=',tign10,&
! 'tign01=',tign01,'tign11=',tign11
!end if
!*** case0 Do nothing
if ( lfn00>=0 .and. lfn10>=0 .and. lfn01>=0 .and. lfn11>=0 ) then
out = 1.0 ! fuel_left, no burning
!*** case4 all four coners are burning
else if (lfn00<=0 .and. lfn10<=0 .and. lfn01<=0 .and. lfn11<=0) then
! least squares, A matrix for points
mA(1,1)=0.0
mA(2,1)=1.0
mA(3,1)=0.0
mA(4,1)=1.0
mA(1,2)=0.0
mA(2,2)=0.0
mA(3,2)=1.0
mA(4,2)=1.0
mA(1,3)=1.0
mA(2,3)=1.0
mA(3,3)=1.0
mA(4,3)=1.0
! D vector, time from ignition
vecD(1)=time_now-tign00
vecD(2)=time_now-tign10
vecD(3)=time_now-tign01
vecD(4)=time_now-tign11
! B matrix, weights
do kk=1,4
do ll=1,4
mB(kk,ll)=0.0
end do
mB(kk,kk)=2.0
end do
! set the m,n,p
n=4 ! rows of matrix A and B
m=3 ! columns of matrix A
p=4 ! columns of matrix B
! call least squqres in LAPACK
call SGGGLM(N,M,P,mA,LDA,mB,LDB,vecD,vecX,vecY, &
WORK,LWORK,INFO)
rnorm=snrm2(p,vecY,1)
! integrate
u(1)=-vecX(1)/fuel_time_cell
u(2)=-vecX(2)/fuel_time_cell
u(3)=-vecX(3)/fuel_time_cell
!fuel_burn(i,j)=1-exp(u(3))*intexp(u(1)*dx)*intexp(u(2)*dy)
s1=u(1)
s2=u(2)
out=1-exp(u(3))*intexp(s1)*intexp(s2)
!print *,'intexp
if ( out<0 .or. out>1.0 ) then
print *,'case4, out should be between 0 and 1'
end if
!*** case 1,2,3
else
! set xx, yy for the coner points
! move these values out of i and j loop to speed up
xx(1) = -0.5
xx(2) = 0.5
xx(3) = 0.5
xx(4) = -0.5
xx(5) = -0.5
yy(1) = -0.5
yy(2) = -0.5
yy(3) = 0.5
yy(4) = 0.5
yy(5) = -0.5
lfn(1)=lfn00
lfn(2)=lfn10
lfn(3)=lfn11
lfn(4)=lfn01
lfn(5)=lfn00
tign(1)=t00
tign(2)=t10
tign(3)=t11
tign(4)=t01
tign(5)=t00
npoint = 0 ! number of points in polygon
!print *,'time_now=',time_now
!print *,'lfn00=',lfn00,'lfn10=',lfn10,&
! 'lfn01=',lfn01,'lfn11=',lfn11
!print *,'t00=',t00,'t10=',t10,&
! 't01=',t01,'t11=',t11
do k=1,4
lfn0=lfn(k )
lfn1=lfn(k+1)
if ( lfn0 <= 0.0 ) then
npoint = npoint + 1
xylist(npoint,1)=xx(k)
xylist(npoint,2)=yy(k)
tlist(npoint)=-tign(k)
llist(npoint)=lfn0
end if
if ( lfn0*lfn1 < 0 ) then
npoint = npoint + 1
tt=lfn0/(lfn0-lfn1)
x0=xx(k)+( xx(k+1)-xx(k) )*tt
y0=yy(k)+( yy(k+1)-yy(k) )*tt
xylist(npoint,1)=x0
xylist(npoint,2)=y0
tlist(npoint)=0 ! on fireline
llist(npoint)=0
end if
end do
! make the list circular
tlist(npoint+1)=tlist(1)
llist(npoint+1)=llist(1)
xylist(npoint+1,1)=xylist(1,1)
xylist(npoint+1,2)=xylist(1,2)
!* least squares, A matrix for points
do kk=1,npoint
mA(kk,1)=xylist(kk,1)
mA(kk,2)=xylist(kk,2)
mA(kk,3)=1.0
vecD(kk)=tlist(kk) ! D vector,time from ignition
end do
! B matrix, weights
do kk=1,ldb
do ll=1,pmax
mB(kk,ll)=0.0 ! clear
end do
end do
do kk=1,npoint
mb(kk,kk) = 10 ! large enough
do ll=1,npoint
if ( kk .ne. ll ) then
dist = sqrt( (xylist(kk,1)-xylist(ll,1))**2+ &
(xylist(kk,2)-xylist(ll,2))**2 )
mB(kk,kk)=min( mB(kk,kk) , dist )
end if
end do !ll
mB(kk,kk)=mB(kk,kk)+1.
end do ! kk
! set the m,n,p
n=npoint ! rows of matrix A and B
m=3 ! columns of matrix A
p=npoint ! columns of matrix B
!* call least squqres in LAPACK
call SGGGLM(N,M,P,mA,LDA,mB,LDB,vecD,vecX,vecY, &
WORK,LWORK,INFO)
!print *,'after LS in case3'
!print *,'vecX from LS',vecX
!print *,'tign inputed',tign00,tign10,tign11,tign01
rnorm=snrm2(p,vecY,1)
u(1)=vecX(1)
u(2)=vecX(2)
u(3)=vecX(3)
! rotate to gradient on x only
nr = sqrt(u(1)**2+u(2)**2)
c = u(1)/nr
s = u(2)/nr
mQ(1,1)=c
mQ(1,2)=s
mQ(2,1)=-s
mQ(2,2)=c
! mat vec multiplication
call matvec(mQ,2,2,u,3,ut,2,2,2)
errQ = ut(2) ! should be zero
ae = -ut(1)/fuel_time_cell
ce = -u(3)/fuel_time_cell
cet=ce!keep ce
call matmatp(xylist,8,2,mQ,2,2,xytlist,8,2,npoint+1,2,2)
call sortxt( xytlist, 8,2, xt,8,npoint )
out=0.0
aupp=0.0
cupp=0.0
alow=0.0
clow=0.0
do k=1,npoint-1
xt1=xt(k)
xt2=xt(k+1)
upper=0
lower=0
ah=0
ch=0
if ( xt2-xt1 > eps*100 ) then
do ss=1,npoint
xts=xytlist(ss,1)
yts=xytlist(ss,2)
xte=xytlist(ss+1,1)
yte=xytlist(ss+1,2)
if ( (xts>xt1 .and. xte>xt1) .or. &
(xts<xt2 .and. xte<xt2) ) then
aa = 0 ! do nothing
cc = 0
else
aa = (yts-yte)/(xts-xte)
cc = (xts*yte-xte*yts)/(xts-xte)
if (xte<xts) then
aupp = aa
cupp = cc
ah=ah+aa
ch=ch+cc
upper=upper+1
else
alow = aa
clow = cc
lower=lower+1
end if
end if!(xts>xt1 .and. xte>xt1)
end do ! ss
ce=cet !use stored ce
if (ae*xt1+ce > 0 ) then
ce=ce-(ae*xt1+ce)!shift small amounts exp(-**)
end if
if (ae*xt2+ce > 0) then
ce=ce-(ae*xt2+ce)
end if
ah = aupp-alow
ch = cupp-clow
! integrate (ah*x+ch)*(1-exp(ae*x+ce) from xt1 to xt2
! numerically sound for ae->0, ae -> infty
! this can be important for different model scales
! esp. if someone runs the model in single precision!!
! s1=int((ah*x+ch),x,xt1,xt2)
s1 = (xt2-xt1)*((1./2.)*ah*(xt2+xt1)+ch)
! s2=int((ch)*(-exp(ae*x+ce)),x,xt1,xt2)
ceae=ce/ae;
s2 = -ch*exp(ae*(xt1+ceae))*(xt2-xt1)*intexp(ae*(xt2-xt1))
! s3=int((ah*x)*(-exp(ae*x+ce)),x,xt1,xt2)
! s3=int((ah*x)*(-exp(ae*(x+ceae))),x,xt1,xt2)
! expand in Taylor series around ae=0
! collect(expand(taylor(int(x*(-exp(ae*(x+ceae))),x,xt1,xt2)*ae^2,ae,4)/ae^2),ae)
! =(1/8*xt1^4+1/3*xt1^3*ceae+1/4*xt1^2*ceae^2-1/8*xt2^4-1/3*xt2^3*ceae-1/4*xt2^2*ceae^2)*ae^2
! + (-1/3*xt2^3-1/2*xt2^2*ceae+1/3*xt1^3+1/2*xt1^2*ceae)*ae
! + 1/2*xt1^2-1/2*xt2^2
!
! coefficient at ae^2 in the expansion, after some algebra
a2=(xt1-xt2)*((1./4.)*(xt1+xt2)*ceae**2+(1./3.)* &
(xt1**2+xt1*xt2+xt2**2)*ceae+(1./8.)* &
(xt1**3+xt1*(xt2**2)+xt1**2*xt2+xt2**3))
d=(ae**4)*a2
if (abs(d)>eps) then
! since ae*xt1+ce<=0 ae*xt2+ce<=0 all fine for large ae
! for ae, ce -> 0 rounding error approx eps/ae^2
s3=( exp(ae*(xt1+ceae))*(ae*xt1-1)-&
exp(ae*(xt2+ceae))*(ae*xt2-1) )/(ae**2)
!we do not worry about rounding as xt1 -> xt2, then s3 -> 0
else
! coefficient at ae^1 in the expansion
a1=(xt1-xt2)*((1./2.)*ceae*(xt1+xt2)+(1./3.)*&
(xt1**2+xt1*xt2+xt2**2))
! coefficient at ae^0 in the expansion for ae->0
a0=(1./2.)*(xt1-xt2)*(xt1+xt2)
s3=a0+a1*ae+a2*ae**2; ! approximate the integral
end if
s3=ah*s3
out=out+s1+s2+s3
out=1-out !fuel_left
if(out<0 .or. out>1) then
print *,':fuel_fraction should be between 0 and 1'
!print *, 'eps= ', eps
!print *, 'xt1= ', xt1, 'xt2= ', xt2
!print *,'ae= ',ae,'ce= ',ce,'ah= ',ah,'ch= ',ch
!print *,'a0= ', a0,'a1= ', a1,'a2= ', a2
print *,'s1= ', s1,'s2= ', s2,'s3= ', s3
print *,':fuel_fraction =',out
end if!print
end if
end do ! k
end if ! if case0, elseif case4 ,else case123
fuel_left_cell = out
end function fuel_left_cell
!
!****************************************
!
real function intexp(ab)
implicit none
real::ab
!calls
intrinsic epsilon
real, parameter:: zero=0.,one=1.,eps=epsilon(zero)
!eps = 2.2204*(10.0**(-8))!from matlab
if ( eps < abs(ab)**3/6. ) then
intexp=(exp(ab)-1)/ab
else
intexp=1+ab/2.
end if
end function
!
!****************************************
!
subroutine sortxt(xytlist,nrow,ncolumn,xt,nxt,nvec)
implicit none
integer::nrow,ncolumn,nxt,nvec
real,dimension(nrow,ncolumn)::xytlist
real,dimension(nxt)::xt
integer::i,j
real::temp
do i=1,nvec
xt(i)=xytlist(i,1)
end do
do i=1,nvec-1
do j=i+1,nvec
if ( xt(i) > xt(j) ) then
temp = xt(i)
xt(i)=xt(j)
xt(j)=temp
end if
end do
end do
end subroutine !sortxt
!
!****************************************
!
subroutine matvec(A,m,n,V,nv,out,nout,nrow,ncolumn)
implicit none
integer::m,n,nv,nout,nrow,ncolumn
real,dimension(m,n)::A ! allocated m by n
real,dimension(nv)::V ! allocated nv
real,dimension(nout)::out! allocated nout
integer::i,j
do i=1,nrow
out(i)=0.0
do j=1,ncolumn
out(i)=out(i)+A(i,j)*V(j)
end do
end do
end subroutine
!
!****************************************
!
subroutine matmatp(A,mA,nA,B,mB,nB,C,mC,nC,nrow,ncolumn,nP)
implicit none
integer::mA,nA,mB,nB,mC,nC,nrow,ncolumn,nP
real,dimension(mA,nA)::A ! allocated m by n
real,dimension(mB,nB)::B ! allocated m by n
real,dimension(mC,nC)::C ! allocated m by n
integer::i,j,k
do i=1,nrow
do j=1,ncolumn
C(i,j)=0.0
do k=1,nP
C(i,j)=C(i,j)+A(i,k)*B(j,k) ! B'
end do
end do
end do
end subroutine
!
!****************************************
!
subroutine prop_ls( id, & ! for debug
ids,ide,jds,jde, & ! domain dims
ims,ime,jms,jme, & ! memory dims
its,ite,jts,jte, & ! tile dims
ts,dt,dx,dy, & ! scalars in
tbound, & ! scalars out
lfn_in,lfn_out,tign & ! arrays inout
#ifdef SPEED_VARS_ARGS /* extra arguments for speed_func */
#include SPEED_VARS_ARGS
#else
#include "fr_sfire_params_args.h"
#endif
)
implicit none
!*** purpose: advance level function in time
! Jan Mandel August 2007 - February 2008
!*** description
!
! Propagation of closed curve by a level function method. The level function
! lfn is defined by its values at the nodes of a rectangular grid.
! The area where lfn < 0 is inside the curve. The curve is
! described implicitly by lfn=0. Points where the curve intersects gridlines
! can be found by linear interpolation from nodes.
!
! The level function is advanced from time ts to time ts + dt.
!
! The level function should be initialized to (an approximation of) the signed
! distance from the curve. If the initial curve is a circle, the initial level
! function is simply the distance from the center minus the radius.
!
! The curve moves outside with speed given by function speed_func.
!
! Method: Godunov method for the normal motion. The timestep is checked for
! CFL condition. For a straight segment in a constant field and locally linear
! level function, the method reduces to the exact normal motion. The advantage of
! the level set method is that it treats automatically special cases such as
! the curve approaching itself and merging components of the area inside the curve.
!
! Based on S. Osher and R. Fedkiw, Level set methods and dynamic implicit surfaces,
! Springer, 2003, Sec. 6.4, as implemented in toolboxLS for Matlab by
! I. Mitchell, A toolbox of Level Set Methods (Version 1.1), TR-2007-11,
! Dept. Computer Science, University of British Columbia, 2007
! http://www.cs.ubc.ca/\~mitchell/Toolbo\LS
!
!*** arguments
! id in unique identification for prints and dumps
! ids,ide,jds,jde in domain dimensions
! ims,ime,jms,jme in memory dimensions
! its,ite,jts,jte in tile dimensions
! ts in start time
! dt in time step
! dx,dy in grid spacing
! lfn_in,lfn_out inout,out the level set function at nodes
! tign inout the ignition time at nodes
! The dimensions are cell-based, the nodal value is associated with the south-west corner.
! The whole computation is on domain indices ids:ide+1,jds:jde+1.
!
! The region where new lfn and tign are computed is the tile its:ite,jts:jte
! except when the tile is at domain upper boundary, an extra band of points is added:
! if ite=ide then region goes up to ite+1, if jte=jde then region goes up to jte+1.
! The time step requires values from 2 rows of nodes beyond the region except when at the
! domain boundary one-sided derivatives are used. This is implemented by extending the input
! beyond the domain boundary so sufficient memory bounds must be allocated.
! The update on all tiles can be done in parallel. To avoid the race condition (different regions
! of the same array updated by different threads), the in and out versions of the
! arrays lft and tign are distinct. If the time step dt is larger
! that the returned tbound, the routine should be called again with timestep td<=tbound, and then
! having distinct inputs and outputs comes handy.
!*** calls
!
! tend_ls
!
integer,intent(in)::id,ims,ime,jms,jme,ids,ide,jds,jde,its,ite,jts,jte
real,dimension(ims:ime,jms:jme),intent(inout)::lfn_in,tign
real,dimension(ims:ime,jms:jme),intent(out)::lfn_out
real,intent(in)::dx,dy,ts,dt
real,intent(out)::tbound
#ifdef SPEED_VARS_DECL /* extra arguments for speed_func */
#include SPEED_VARS_DECL
#else
#include "fr_sfire_params_decl.h"
#endif
!*** local
! arrays
#define IMTS its-1
#define IMTE ite+1
#define JMTS jts-1
#define JMTE jte+1
real,dimension(IMTS:IMTE,JMTS:JMTE):: tend, lfn1 ! region-sized with halo
! scalars
real::grad2,rr,tbound2
real, parameter::a=0.5, a1=1.0-a ! a=0 euler, a=0.5 heun
real::gradx,grady,aspeed,err,aerr,time_now
integer::ihs,ihe,jhs,jhe
integer::i,j,its1,ite1,jts1,jte1
character(len=128)msg
integer::nfirenodes,nfireline
real::sum_err,min_err,max_err,sum_aerr,min_aerr,max_aerr
! constants
integer,parameter :: mstep=1000, printl=1
real, parameter:: zero=0.,one=1.,eps=epsilon(zero),tol=100*eps, &
safe=2.,rmin=safe*tiny(zero),rmax=huge(zero)/safe
! f90 intrinsic function
intrinsic max,min,sqrt,nint,epsilon,tiny,huge
!*** executable
write(msg,'(5(a,i5))')'prop_ls:',id,' tile ',its,':',ite,',',jts,':',jte
call message(msg)
! tend = F(lfn)
ihs=max(its-1,ids) ! compute tend one beyond the tile but not outside the domain
ihe=min(ite+1,ide)
jhs=max(jts-1,jds)
jhe=min(jte+1,jde)
#ifdef DEBUG_OUT
call write_array_m(ihs,ihe,jhs,jhe,ims,ime,jms,jme,lfn_in,'lfn_in',id)
#endif
! NOTE: tend_ls will do relection to one node strip at domain boundaries
! so that it can compute gradient at domain boundaries.
! To avoid copying, lfn_in is declared inout.
! At tile boundaries that are not domain boundaries values of lfn_in two nodes
! outside of the tile are needed.
call tend_ls(&
ims,ime,jms,jme, & ! memory dims for lfn_in
IMTS,IMTE,JMTS,JMTE, & ! memory dims for tend
ids,ide,jds,jde, & ! domain dims - where lfn exists
ihs,ihe,jhs,jhe, & ! where tend computed
ts,dx,dy, & ! scalars in
lfn_in, & ! arrays in
tbound, & ! scalars out
tend & ! arrays out
#ifdef SPEED_VARS_ARGS /* extra arguments for speed_func */
#include SPEED_VARS_ARGS
#else
#include "fr_sfire_params_args.h"
#endif
)
#ifdef DEBUG_OUT
call write_array_m(ihs,ihe,jhs,jhe,IMTS,IMTE,JMTS,JMTE,tend,'tend1',id)
#endif
! Euler method, the half-step, same region as ted
do j=jhs,jhe
do i=ihs,ihe
lfn1(i,j) = lfn_in(i,j) + dt*tend(i,j)
enddo
enddo
! tend = F(lfn1) on the tile (not beyond)
call tend_ls( &
IMTS,IMTE,JMTS,JMTE, & ! memory dims for lfn
IMTS,IMTE,JMTS,JMTE, & ! memory dims for tend
ids,ide,jds,jde, & ! domain dims - where lfn exists
its,ite,jts,jte, & ! tile dims - where is tend computed
ts+dt,dx,dy, & ! scalars in
lfn1, & ! arrays in
tbound2, & ! scalars out
tend & ! arrays out
#ifdef SPEED_VARS_ARGS /* extra arguments for speed_func */
#include SPEED_VARS_ARGS
#else
#include "fr_sfire_params_args.h"
#endif
)
call print_2d_stats(its,ite,jts,jte,IMTS,IMTE,JMTS,JMTE,tend,'prop_ls:lvl fcn tend(1/s)')
tbound=min(tbound,tbound2)
write(msg,'(a,f10.2,4(a,f7.2))')'prop_ls: time',ts,' dt=',dt,' bound',min(tbound,999.99), &
' dx=',dx,' dy=',dy
call message(msg)
if(dt>tbound)then
write(msg,'(2(a,f10.2))')'prop_ls: WARNING: time step ',dt, &
' > bound =',tbound
call message(msg)
endif
! combine lfn1 and lfn_in + dt*tend -> lfn_out
do j=jts,jte
do i=its,ite
lfn_out(i,j) = a1*lfn1(i,j) + a*(lfn_in(i,j) + dt*tend(i,j))
enddo
enddo
! compute ignition time by interpolation
! the node was not burning at start but it is burning at end
! interpolate from the level functions at start and at end
! lfn_in is the level set function value at time ts
! lfn_out is the level set function value at time ts+dt
! 0 is the level set function value at time tign(i,j)
! thus assuming the level function is approximately linear =>
! tign(i,j)= ts + ((ts + td) - ts) * lfn_in / (lfn_in - lfn_out)
! = ts + dt * lfn_in / (lfn_in - lfn_out)
time_now=ts+dt