debug
[wrf-fire-matlab.git] / femwind / fortran / module_boundary_conditions.f90
blob1c0c5311f5aa4693c1cba34721732e0451163724
1 module module_boundary_conditions
3 use module_utils
5 contains
7 subroutine ndt_boundary_conditions( &
8 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
9 ifms, ifme, kfms, kfme, jfms, jfme, &
10 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
11 ifts, ifte, kfts, kfte, jfts, jfte, &
12 kmat)
14 implicit none
16 !*** arguments
18 integer, intent(in):: &
19 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
20 ifms, ifme, kfms, kfme, jfms, jfme, &
21 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
22 ifts, ifte, kfts, kfte, jfts, jfte
23 integer, parameter:: msize = 14
24 real, intent(inout), dimension(ifms:ifme,kfms:kfme,jfms:jfme,msize):: kmat ! global stiffness matrix
26 !*** local
28 real:: s
29 integer:: i,j,k,ie,je,ke
31 !** executable
33 ! loop upper bounds
34 ie = snode(ifte,ifde,+1)
35 je = snode(jfte,jfde,+1)
36 ke = snode(kfte,kfde,+1)
38 ! scale
39 s=0.
41 do j=jfts,je
42 do k=kfts,ke
43 do i=ifts,ie
44 s=max(s,abs(kmat(i ,k ,j , 1)))
45 enddo
46 enddo
47 enddo
49 do j=jfts,je
50 do k=kfts,ke
51 do i=ifts,ie
52 ! not efficient but will be executed only once
53 if(i.eq.ifds.or.i.eq.ifde+1.or.j.eq.jfds.or.j.eq.jfde+1.or.k.eq.kfde+1)then
54 ! replace row/col (i,k,j) by scaled identity
55 print *,'replacing row/col ',i,j,k,' by ',s, '* identity'
56 kmat(i-1,k-1,j-1,14)=0.
57 kmat(i ,k-1,j-1,13)=0.
58 kmat(i+1,k-1,j-1,12)=0.
59 kmat(i-1,k-1,j ,11)=0.
60 kmat(i ,k-1,j ,10)=0.
61 kmat(i+1,k-1,j , 9)=0.
62 kmat(i-1,k-1,j+1, 8)=0.
63 kmat(i ,k-1,j+1, 7)=0.
64 kmat(i+1,k-1,j+1, 6)=0.
65 kmat(i-1,k ,j-1, 5)=0.
66 kmat(i ,k ,j-1, 4)=0.
67 kmat(i+1,k ,j-1, 3)=0.
68 kmat(i-1,k ,j , 2)=0.
69 kmat(i ,k ,j , 1)=s
70 kmat(i ,k ,j , 2)=0.
71 kmat(i ,k ,j , 3)=0.
72 kmat(i ,k ,j , 4)=0.
73 kmat(i ,k ,j , 5)=0.
74 kmat(i ,k ,j , 6)=0.
75 kmat(i ,k ,j , 7)=0.
76 kmat(i ,k ,j , 8)=0.
77 kmat(i ,k ,j , 9)=0.
78 kmat(i ,k ,j ,10)=0.
79 kmat(i ,k ,j ,11)=0.
80 kmat(i ,k ,j ,12)=0.
81 kmat(i ,k ,j ,13)=0.
82 kmat(i ,k ,j ,14)=0.
83 endif
84 enddo
85 enddo
86 enddo
88 end subroutine ndt_boundary_conditions
90 subroutine vec_boundary_conditions( &
91 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
92 ifms, ifme, kfms, kfme, jfms, jfme, &
93 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
94 ifts, ifte, kfts, kfte, jfts, jfte, &
97 implicit none
99 !*** arguments
101 integer, intent(in):: &
102 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
103 ifms, ifme, kfms, kfme, jfms, jfme, &
104 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
105 ifts, ifte, kfts, kfte, jfts, jfte
106 integer, parameter:: msize = 14
107 real, intent(inout), dimension(ifms:ifme,kfms:kfme,jfms:jfme):: F ! corner-based scalar field
109 !*** local
111 integer:: i,j,k,ie,je,ke
113 !** executable
115 ! loop upper bounds
116 ie = snode(ifte,ifde,+1)
117 je = snode(jfte,jfde,+1)
118 ke = snode(kfte,kfde,+1)
120 !call write_array(F(ifts:ie,kfts:ke,jfts:je),'F_bc_in')
122 do j=jfts,je
123 do k=kfts,ke
124 do i=ifts,ie
125 ! not efficient, change later
126 if(i.eq.ifds.or.i.eq.ifde+1.or.j.eq.jfds.or.j.eq.jfde+1.or.k.eq.kfde+1)then
127 F(i,k,j)=0.
128 endif
129 enddo
130 enddo
131 enddo
133 !call write_array(F(ifts:ie,kfts:ke,jfts:je),'F_bc_out')
135 end subroutine vec_boundary_conditions
136 end module module_boundary_conditions