femwind_rate_test and femwind_run_fortran_rate_test both exiting in iteration 20...
[wrf-fire-matlab.git] / femwind / fortran / module_f_assembly.f90
blobd575c3d5569d0471f5b4b6770dd3198af2a4f1e1
1 module module_f_assembly
2 use module_hexa
3 use module_utils
5 contains
6 subroutine f_assembly( &
7 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
8 ifms, ifme, kfms, kfme, jfms, jfme, &
9 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
10 ifts, ifte, kfts, kfte, jfts,jfte, &
11 A, X, Y, Z, Xu0, Yu0, Zu0, & !Input from femwind, U0, V0, W0 not used in hexa to construct Kloc or JG
12 F) !Global load vector output
14 implicit none
17 !*** arguments
19 integer, intent(in):: &
20 ifds, ifde, kfds,kfde, jfds, jfde, & ! fire grid dimensions
21 ifms, ifme, kfms,kfme, jfms, jfme, &
22 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
23 ifts, ifte, kfts, kfte, jfts,jfte
28 real, intent(in), dimension(3,3):: A
29 real, intent(in), dimension(ifms:ifme, kfms:kfme, jfms:jfme):: X,Y,Z, & !spatial grid at cornersw
30 Xu0,Yu0,Zu0 !wind vector at midpoint
31 !Input for hexa
33 !integer,intent(in)::F_dim
35 real,intent(out), dimension(ifms:ifme, kfms:kfme, jfms:jfme)::F
37 !*** local
39 integer:: ie1, ie2, ie3, ic1, ic2, ic3, iloc, k1, k2, k3, id1, id2, id3
40 real :: Kloc(8,8), Floc(8), Jg(8,3)
41 real :: Xloc(3,8), u0loc(3)
42 integer :: kglo(8)
43 !*** u0loc is an input for module_hexa, but is not used to construct K. Do I need to define this?
44 !*** integer, dimension(3,1,1), save ::iflags = reshape((/1,0,1/),(/3,1,1/)) !define iflags to construct JG and Kloc in hexa
46 call write_array(reshape(A,(/3,3,1/)),'A_in')
47 call write_array(X(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'X_in')
48 call write_array(Y(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'Y_in')
49 call write_array(Z(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'Z_in')
50 call write_array(Xu0(ifts: ifte, kfts: kfte, jfts:jfte),'u0_in')
51 call write_array(Yu0(ifts: ifte, kfts: kfte, jfts:jfte),'v0_in')
52 call write_array(Zu0(ifts: ifte, kfts: kfte, jfts:jfte),'w0_in')
54 Xloc = 9999.
55 u0loc = 0.
56 F = 0.
58 !** executable
59 do ie2=jfts,jfte
60 do ie3=kfts,kfte
61 do ie1=ifts,ifte
62 do ic2=0,1
63 do ic3=0,1
64 do ic1=0,1
65 iloc=1+ic1+2*(ic2+2*ic3) !local index of the node in the element
66 Xloc(1,iloc)=X(ie1 + ic1, ie3 + ic3, ie2 + ic2)
67 Xloc(2,iloc)=Y(ie1 + ic1, ie3 + ic3, ie2 + ic2)
68 Xloc(3,iloc)=Z(ie1 + ic1, ie3 + ic3, ie2 + ic2)
69 enddo
70 enddo
71 enddo
72 u0loc(1) = Xu0(ie1,ie3,ie2)
73 u0loc(2) = Yu0(ie1,ie3,ie2)
74 u0loc(3) = Zu0(ie1,ie3,ie2)
76 call hexa(A,Xloc,u0loc,Kloc,Floc,Jg,2)
77 do id2 = 0,1
78 do id3 = 0,1
79 do id1 = 0,1
80 iloc=1+id1+2*(id2+2*id3) !local index of the node in the element
81 k1 = ie1+id1
82 k2 = ie2+id2
83 k3 = ie3+id3
84 F(k1,k3,k2) = F(k1,k3,k2) + Floc(iloc)
85 enddo
86 enddo
87 end do
88 enddo
89 enddo
90 enddo
92 call write_array(F(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'F_out')
94 end subroutine f_assembly
95 end module module_f_assembly