Merge branch 'fixf'
[wrf-fire-matlab.git] / femwind / fortran / module_femwind.f90
bloba151b15fcc6fd38a6a63dcb02e257f3c17c4a9d0
1 ! **** HOW TO TEST ****
2 ! run femwind_fortran_rate_test
4 ! details:
5 ! called from femwind_test.f90 compiled as femwind_test.exe
6 ! femwind_test.exe is executed from femwind_fortran.m if params.run_fortran
7 ! femwind_fortran.m is called from femwind_main.m
8 ! run femwind_main from e.g. femwind_rate_test or femwind_wrfout_test
9 ! femwind_fortran_rate_test sets the params and calls femwind_rate_test
12 module module_femwind
13 use module_multigrid
14 use module_common
15 use module_utils
16 use module_f_assembly
17 use module_ndt_assembly
18 use module_w_assembly
19 ! use module_ndt_mult
20 ! use module_coarsening
21 use module_boundary_conditions
22 ! use module_sweeps
24 contains
26 subroutine femwind_setup(mg)
27 implicit none
28 ! set up the mg_level structure
29 ! input: mg(1)%X,Y,Z,dx,dy,dz already set
30 !*** arguments
31 type(mg_type),intent(inout)::mg(:) ! multigrid level
32 !*** local
34 !*** executable
36 call multigrid_setup(mg)
38 end subroutine femwind_setup
40 subroutine femwind_solve(mg, &
41 ifds, ifde, kfds,kfde, jfds, jfde, & ! fire grid dimensions
42 ifms, ifme, kfms,kfme, jfms, jfme, &
43 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
44 ifts, ifte, kfts, kfte, jfts,jfte, &
45 u0, v0, w0, & ! input
46 u, v, w, rate) ! output
48 implicit none
50 !*** arguments
52 type(mg_type),intent(inout)::mg(:) ! multigrid levels
54 integer, intent(in):: &
55 ifds, ifde, kfds,kfde, jfds, jfde, & ! fire grid dimensions
56 ifms, ifme, kfms,kfme, jfms, jfme, &
57 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
58 ifts, ifte, kfts, kfte, jfts,jfte
60 real, intent(in), dimension(ifms:ifme, kfms:kfme, jfms:jfme):: u0, v0, w0 !initial wind vector at midpoints
61 real,intent(out), dimension(ifms:ifme, kfms:kfme, jfms:jfme):: u, v, w !mass consistent wind at midpoints
62 real,intent(out):: rate
64 !*** local
66 !*** executable
68 ! f = div(u0)
69 ! f = f_assembly_fortran(A,X,U0,lambda,params);
71 print *,'femwind solve start'
72 print *,'calling f_assembly'
73 call f_assembly( &
74 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
75 ifms, ifme, kfms, kfme, jfms, jfme, &
76 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
77 ifts, ifte, kfts, kfte, jfts, jfte, &
78 params%A, mg(1)%X, mg(1)%Y, mg(1)%Z, &
79 u0, v0, w0, &
80 mg(1)%f) !Global load vector output
82 print *,'calling vec_boundary_conditions'
83 call vec_boundary_conditions( &
84 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
85 ifms, ifme, kfms, kfme, jfms, jfme, &
86 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
87 ifts, ifte, kfts, kfte, jfts, jfte, &
88 mg(1)%f)
90 print *,'calling multigrid_cycle'
91 call multigrid_cycle(mg,1,rate) ! start multigrid from level 1
92 print *,'end multigrid_cycle, rate=',rate
94 if(params%debug_level >=0)call write_array(mg(1)%lambda(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'lambda_sol')
96 if(params%debug_level >=0)call write_array(u0(ifts: ifte, kfts: kfte, jfts:jfte),'u0_sol')
97 if(params%debug_level >=0)call write_array(v0(ifts: ifte, kfts: kfte, jfts:jfte),'v0_sol')
98 if(params%debug_level >=0)call write_array(w0(ifts: ifte, kfts: kfte, jfts:jfte),'w0_sol')
99 if(params%debug_level >=0)call write_array(mg(1)%X(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'X_sol')
100 if(params%debug_level >=0)call write_array(mg(1)%Y(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'Y_sol')
101 if(params%debug_level >=0)call write_array(mg(1)%Z(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'Z_sol')
103 print *,'calling w_assembly'
104 call w_assembly( &
105 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
106 ifms, ifme, kfms, kfme, jfms, jfme, &
107 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
108 ifts, ifte, kfts, kfte, jfts, jfte, &
109 mg(1)%lambda, u0, v0, w0, &
110 params%A, mg(1)%X, mg(1)%Y, mg(1)%Z, & !Input from femwind, u0, v0, w0, Spatial Grid Data
111 u, v, w) ! final output
113 print *,'end femwind_solve'
115 if(params%debug_level >=0)call write_array(u(ifts: ifte, kfts: kfte, jfts:jfte),'u_sol')
116 if(params%debug_level >=0)call write_array(v(ifts: ifte, kfts: kfte, jfts:jfte),'v_sol')
117 if(params%debug_level >=0)call write_array(w(ifts: ifte, kfts: kfte, jfts:jfte),'w_sol')
119 end subroutine femwind_solve
122 end module module_femwind