max_hfx to diagnose computation of fgrnhfx around max
[wrf-fire-matlab.git] / femwind / fortran / module_femwind.f90
blob50c648d8743fc985e624423ce73e55469e6734bf
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 call f_assembly( &
72 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
73 ifms, ifme, kfms, kfme, jfms, jfme, &
74 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
75 ifts, ifte, kfts, kfte, jfts, jfte, &
76 params%A, mg(1)%X, mg(1)%Y, mg(1)%Z, &
77 u0, v0, w0, &
78 mg(1)%f) !Global load vector output
80 call vec_boundary_conditions( &
81 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
82 ifms, ifme, kfms, kfme, jfms, jfme, &
83 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
84 ifts, ifte, kfts, kfte, jfts, jfte, &
85 mg(1)%f)
87 call multigrid_cycle(mg,1,rate) ! start multigrid from level 1
89 if(params%debug_level >=0)call write_array(mg(1)%lambda(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'lambda_sol')
91 if(params%debug_level >=0)call write_array(u0(ifts: ifte, kfts: kfte, jfts:jfte),'u0_sol')
92 if(params%debug_level >=0)call write_array(v0(ifts: ifte, kfts: kfte, jfts:jfte),'v0_sol')
93 if(params%debug_level >=0)call write_array(w0(ifts: ifte, kfts: kfte, jfts:jfte),'w0_sol')
94 if(params%debug_level >=0)call write_array(mg(1)%X(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'X_sol')
95 if(params%debug_level >=0)call write_array(mg(1)%Y(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'Y_sol')
96 if(params%debug_level >=0)call write_array(mg(1)%Z(ifts: ifte+1, kfts: kfte+1, jfts:jfte+1),'Z_sol')
98 call w_assembly( &
99 ifds, ifde, kfds, kfde, jfds, jfde, & ! fire grid dimensions
100 ifms, ifme, kfms, kfme, jfms, jfme, &
101 ifps, ifpe, kfps, kfpe, jfps, jfpe, & ! fire patch bounds
102 ifts, ifte, kfts, kfte, jfts, jfte, &
103 mg(1)%lambda, u0, v0, w0, &
104 params%A, mg(1)%X, mg(1)%Y, mg(1)%Z, & !Input from femwind, u0, v0, w0, Spatial Grid Data
105 u, v, w) ! final output
107 if(params%debug_level >=0)call write_array(u(ifts: ifte, kfts: kfte, jfts:jfte),'u_sol')
108 if(params%debug_level >=0)call write_array(v(ifts: ifte, kfts: kfte, jfts:jfte),'v_sol')
109 if(params%debug_level >=0)call write_array(w(ifts: ifte, kfts: kfte, jfts:jfte),'w_sol')
111 end subroutine femwind_solve
114 end module module_femwind