ndt_boundary_conditions_test compiles
[wrf-fire-matlab.git] / femwind / fortran / module_boundary_conditions.f90
blobd63ff0a219c2d8578ca3c8a402e178f0b0c0bed5
1 module module_boundary_conditions
3 contains
5 subroutine ndt_boundary_conditions( &
6 ifds, ifde, kfds,kfde, jfds, jfde, & ! fire grid dimensions
7 ifms, ifme, kfms,kfme, jfms, jfme, &
8 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
9 ifts, ifte, kfts, kfte, jfts,jfte, &
10 kmat)
12 implicit none
14 !*** arguments
16 integer, intent(in):: &
17 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire domain bounds
18 ifms, ifme, kfms, kfme, jfms, jfme, & ! fire memory bounds
19 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
20 ifts, ifte, kfts, kfte, jfts,jfte ! fire tile bounds
22 integer, parameter:: msize = 14
23 real, intent(inout), dimension(ifms:ifme,kfms:kfme,jfms:jfme,msize):: kmat ! global stiffness matrix
25 !*** local
27 real:: s
28 integer:: i,j,k
30 !** executable
32 ! scale
33 s=0.
35 do j=jfts,jfte
36 do k=kfts,kfte
37 do i=ifts,ifte
38 s=max(s,abs(kmat(i ,k ,j , 1)))
39 enddo
40 enddo
41 enddo
43 do j=jfts,jfte
44 do k=kfts,kfte
45 do i=ifts,ifte
46 ! not efficient but will be executed only once
47 if(i.eq.ifds.or.i.eq.ifde.or.j.eq.jfds.or.j.eq.jfde.or.k.eq.kfde)then
48 ! replace the row/col (i,k,j) by scaled identity
49 kmat(i-1,k-1,j-1,14)=0.
50 kmat(i ,k-1,j-1,13)=0.
51 kmat(i+1,k-1,j-1,12)=0.
52 kmat(i-1,k-1,j ,11)=0.
53 kmat(i ,k-1,j ,10)=0.
54 kmat(i+1,k-1,j , 9)=0.
55 kmat(i-1,k-1,j+1, 8)=0.
56 kmat(i ,k-1,j+1, 7)=0.
57 kmat(i+1,k-1,j+1, 6)=0.
58 kmat(i-1,k ,j-1, 5)=0.
59 kmat(i ,k ,j-1, 4)=0.
60 kmat(i+1,k ,j-1, 3)=0.
61 kmat(i-1,k ,j , 2)=0.
62 kmat(i ,k ,j , 1)=s
63 kmat(i ,k ,j , 2)=0.
64 kmat(i ,k ,j , 3)=0.
65 kmat(i ,k ,j , 4)=0.
66 kmat(i ,k ,j , 5)=0.
67 kmat(i ,k ,j , 6)=0.
68 kmat(i ,k ,j , 7)=0.
69 kmat(i ,k ,j , 8)=0.
70 kmat(i ,k ,j , 9)=0.
71 kmat(i ,k ,j ,10)=0.
72 kmat(i ,k ,j ,11)=0.
73 kmat(i ,k ,j ,12)=0.
74 kmat(i ,k ,j ,13)=0.
75 kmat(i ,k ,j ,14)=0.
76 endif
77 enddo
78 enddo
79 enddo
81 end subroutine ndt_boundary_conditions
83 end module module_boundary_conditions