MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
tests.f90
Go to the documentation of this file.
1 
9 module tests
10  use foaming_globals_m
11  use constants
12  use in_out, only:set_paths,read_inputs,load_old_results
13  use integration, only:bblpreproc,bblinteg
14  implicit none
15  private
16  public onegrowth,eta_rm,bub_vf,secondgrowth,shooting_method,&
17  shooting_method_test
18 contains
19 !********************************BEGINNING*************************************
25 subroutine onegrowth
26  use phys_prop, only:rb_initialized
27  rb_initialized=.false.
28  call set_paths
29  call read_inputs
30  call bblpreproc
31  call bblinteg
32  rb_initialized=.false.
33 end subroutine onegrowth
34 !***********************************END****************************************
35 
36 
37 !********************************BEGINNING*************************************
42 subroutine secondgrowth
43  use globals
44  use fson
45  use fson_value_m, only: fson_value_get
46  real(dp) :: tt,rr
47  type(fson_value), pointer :: json_data
48  ! prepare to work with files
49  call set_paths
50  ! obtain tend
51  write(*,*) 'loading input file ',trim(inputs)
52  json_data => fson_parse(inputs)
53  call fson_get(json_data, "bubbleGrowth.finalTime", tend)
54  ! feed results to function Rb(t)
55  call load_old_results
56  bub_inx=1
57  ! simulation
58  firstrun=.false.
59  call onegrowth
60 end subroutine secondgrowth
61 !***********************************END****************************************
62 
63 
64 !********************************BEGINNING*************************************
69 real(dp) function nextradius(n,rad_ini)
70  use globals
71  integer, intent(in) :: n
72  real(dp), dimension(n), intent(in) :: rad_ini
74  r0=rad_ini(1)
75  call onegrowth
76  nextradius=radius
77 end function nextradius
78 !***********************************END****************************************
79 
80 
81 !********************************BEGINNING*************************************
86 subroutine rad_residual(n,x,fvec,iflag)
87  use globals
88  integer, intent(in) :: n
89  integer, intent(inout) :: iflag
90  real(dp), dimension(n), intent(in) :: x
91  real(dp), dimension(n), intent(out) :: fvec
92  fvec(1)=sqrt((nextradius(n,x)-goalradius)**2)
93  print*, x(1),fvec(1)
94 end subroutine rad_residual
95 !***********************************END****************************************
96 
97 
98 !********************************BEGINNING*************************************
103 subroutine shooting_method
105  use solve_nonlin
106  use fson
107  use fson_value_m, only: fson_value_get
108  integer, parameter :: n=1
109  integer :: info
110  real (dp) :: rad_ini,frad_ini,tol=1e-8_dp
111  real (dp), dimension(n) :: x,fvec,diag
112  type(fson_value), pointer :: json_data
113  ! prepare to work with files
114  call set_paths
115  ! obtain tend
116  write(*,*) 'loading input file ',trim(inputs)
117  json_data => fson_parse(inputs)
118  call fson_get(json_data, "initialConditions.bubbleRadius", r0)
119  rad_ini=r0
120  goalradius=bub_rad(size(bub_rad(:,1)),bub_inx+1)
121  firstrun=.true.
122  shooting=.true.
123  x(1)=rad_ini
124  call hbrd(rad_residual,n,x,fvec,epsilon(pi),tol,info,diag)
125  if (info /= 1) then
126  write(unit=*, fmt=*) 'Shooting method failed.'
127  write(unit=*, fmt=*) 'Hbrd returned info = ',info
128  stop
129  endif
130  r0=x(1)
131  call onegrowth
132  write(unit=*, fmt=*) 'First guess initial radius = ',rad_ini
133  write(unit=*, fmt=*) 'Converged initial radius = ',r0
134  write(unit=*, fmt=*) 'Goal radius = ',goalradius
135  write(unit=*, fmt=*) 'Final converged radius = ',radius
136 end subroutine shooting_method
137 !***********************************END****************************************
138 
139 
140 !********************************BEGINNING*************************************
145 subroutine shooting_method_test
147  use solve_nonlin
148  integer, parameter :: n=1
149  integer :: info
150  real (dp) :: rad_ini,frad_ini,tol=1e-8_dp
151  real (dp), dimension(n) :: x,fvec,diag
152  ! first simulation
153  firstrun=.true.
154  shooting=.false.
155  call onegrowth
156  rad_ini=r0
157  frad_ini=radius
158  goalradius=radius*(1+1.e-4_dp)
159  shooting=.true.
160  x(1)=rad_ini
161  call hbrd(rad_residual,n,x,fvec,epsilon(pi),tol,info,diag)
162  if (info /= 1) then
163  write(unit=*, fmt=*) 'Shooting method failed.'
164  write(unit=*, fmt=*) 'Hbrd returned info = ',info
165  stop
166  endif
167  r0=x(1)
168  call onegrowth
169  write(unit=*, fmt=*) 'First guess initial radius = ',rad_ini
170  write(unit=*, fmt=*) 'Converged initial radius = ',r0
171  write(unit=*, fmt=*) 'Goal radius = ',goalradius
172  write(unit=*, fmt=*) 'First guess final radius = ',frad_ini
173  write(unit=*, fmt=*) 'Final converged radius = ',radius
174 end subroutine shooting_method_test
175 !***********************************END****************************************
176 
177 
178 
179 !********************************BEGINNING*************************************
183 real(dp) function eta_rm(t)
184  use interpolation
185  real(dp), intent(in) :: t
186  integer :: n
187  integer :: ni=1 !number of points, where we want to interpolate
188  real(dp) :: xi(1) !x-values of points, where we want to interpolate
189  real(dp) :: yi(1) !interpolated y-values
190  xi(1)=t
191  call pwl_interp_1d ( size(etat(:,1)), etat(:,1), etat(:,2), ni, xi, yi )
192  eta_rm=yi(1)
193 endfunction eta_rm
194 !***********************************END****************************************
195 
196 
197 !********************************BEGINNING*************************************
201 real(dp) function bub_vf(t)
202  use interpolation
203  real(dp), intent(in) :: t
204  integer :: ni=1 !number of points, where we want to interpolate
205  real(dp) :: xi(1) !x-values of points, where we want to interpolate
206  real(dp) :: yi(1) !interpolated y-values
207  xi(1)=t
208  call pwl_interp_1d ( size(port(:,1)), port(:,1), port(:,2), ni, xi, yi )
209  bub_vf=yi(1)
210 endfunction bub_vf
211 !***********************************END****************************************
212 end module tests
character(len=99) inputs
input file
Definition: globals.f90:11
real(dp) radius
bubble radius (m)
Definition: globals.f90:49
real(dp) goalradius
final radius we want to achieve
Definition: globals.f90:49
logical shooting
am I using shooting method (t/f)
Definition: globals.f90:22
real(dp) r0
initial radius
Definition: globals.f90:49
namespace with global variables
Definition: globals.f90:8