MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
in_out.f90
Go to the documentation of this file.
1 
7 module in_out
8  implicit none
9  private
10  integer, dimension(:), allocatable :: fi
11  public read_inputs,save_int_header,save_int_step,save_int_close,&
12  load_bubble_growth
13 contains
14 !********************************BEGINNING*************************************
16 subroutine read_inputs
17  use globals
18  use fson
19  use fson_value_m, only: fson_value_get
20  character(len=1024) :: strval
21  type(fson_value), pointer :: json_data
22  json_data => fson_parse("../inputs/wallDrainage_inputs.json")
23  call fson_get(json_data, "initialConditions.centerThickness", hi)
24  call fson_get(json_data, "growthRateModel", growthratemodel)
25  if (growthratemodel=="constantGrowth") then
26  call fson_get(json_data, "growthRate", gr)
27  call fson_get(json_data, "initialConditions.domainSize", rd)
28  call fson_get(json_data, "initialConditions.filmReduction", dstr)
29  call fson_get(json_data, "integration.initialTime", initialtime)
30  call fson_get(json_data, "integration.outerTimeSteps", its)
31  elseif (growthratemodel=="fromFile") then
32  continue
33  else
34  print*, 'unknown growth rate model'
35  stop
36  endif
37  call fson_get(json_data, "physicalProperties.viscosityModel", viscositymodel)
38  if (viscositymodel=="constant") then
39  call fson_get(json_data, "physicalProperties.viscosity", mu)
40  elseif (viscositymodel=="fromFile") then
41  continue
42  else
43  print*, 'unknown viscosity model'
44  stop
45  endif
46  call fson_get(json_data, "physicalProperties.surfaceTension", gam)
47  call fson_get(json_data, "physicalProperties.disjoiningPressure.N", ndp)
48  call fson_get(json_data, "physicalProperties.disjoiningPressure.M", mdp)
49  call fson_get(json_data, "physicalProperties.disjoiningPressure.hst", cdp)
50  call fson_get(json_data, "physicalProperties.disjoiningPressure.C", hdp)
51  call fson_get(json_data, "physicalProperties.disjoiningPressure.B1", bdp)
52  call fson_get(json_data, "integration.timeStep", timestep)
53  call fson_get(json_data, "integration.method", strval)
54  if (strval=="stiff") then
55  int_method=222
56  else
57  print*, "unknown integration method"
58  stop
59  endif
60  call fson_get(json_data, "integration.internalNodes", meshpoints)
61  call fson_get(json_data, "integration.maxInnerTimeSteps", maxts)
62  call fson_get(json_data, "integration.relativeTolerance", int_reltol)
63  call fson_get(json_data, "integration.absoluteTolerance", int_abstol)
64  call fson_get(json_data, "algebraicEquationSolver.tolerance", ae_tol)
65  call fson_get(json_data, "strutFilmParameter", strutfilmparameter)
66  call fson_destroy(json_data)
67 end subroutine read_inputs
68 !***********************************END****************************************
69 
70 
71 !********************************BEGINNING*************************************
73 subroutine save_int_header
74  use ioutils, only: newunit
75  allocate(fi(10))
76  open (unit=newunit(fi(1)), file = 'filmthickness.csv')
77  open (unit=newunit(fi(2)), file = 'results_1d.csv')
78  write(*,'(1x,100a12)') 'time:','dr:','total: ','struts:',&
79  'hmin: ','hloc: ','hcenter: ','havg: '
80  write(unit=fi(2), fmt='(10000a12)') 'time','dr','np','vt','fs',&
81  'hmin','hloc','hcenter','havg'
82 end subroutine save_int_header
83 !***********************************END****************************************
84 
85 
86 !********************************BEGINNING*************************************
88 subroutine save_int_step(y,t)
89  use constants, only: dp
90  use globals, only: dr
91  use phys_prop, only: volume_balance,min_film_thickness
92  real(dp), intent(in) :: t
93  real(dp), dimension(:), intent(in) :: y
94  integer :: neq
95  real(dp) :: fs,vt,hmin,hloc,havg
96  neq=size(y)
97  call volume_balance(y,vt,fs)
98  call min_film_thickness(y,hmin,hloc,havg)
99  write(*,'(100es12.3)') t,dr,vt,fs,hmin,hloc,y(1),havg
100  write(fi(1),"(10000es12.4)") y(1:neq)
101  write(unit=fi(2), fmt='(10000es12.4)') &
102  t,dr,dble(neq),vt,fs,hmin,hloc,y(1),havg
103 end subroutine save_int_step
104 !***********************************END****************************************
105 
106 
107 !********************************BEGINNING*************************************
109 subroutine save_int_close
110  close(fi(1))
111  close(fi(2))
112 end subroutine save_int_close
113 !***********************************END****************************************
114 
115 
116 !********************************BEGINNING*************************************
118 subroutine load_bubble_growth(matrix)
119  use constants, only: dp
120  use ioutils, only: newunit
122  real(dp), dimension(:,:), allocatable, intent(inout) :: matrix
123  integer :: i,j,ios,fi
124  j=0
125  open(newunit(fi),file='../results/bubbleGrowth/bblgr_2_drain.out')
126  do !find number of points
127  read(fi,*,iostat=ios)
128  if (ios/=0) exit
129  j=j+1
130  enddo
131  if (allocated(matrix)) deallocate(matrix)
132  allocate(matrix(j-1,4))
133  rewind(fi)
134  read(fi,*)
135  do i=1,j-1
136  read(fi,*) matrix(i,:)
137  enddo
138  close(fi)
139 end subroutine load_bubble_growth
140 !***********************************END****************************************
141 end module in_out
real(dp) gam
surface tension
Definition: globals.f90:18
integer meshpoints
number of discretization points in space
Definition: globals.f90:13
character(len=1024) growthratemodel
growth rate model (constantGrowth, fromFile)
Definition: globals.f90:10
real(dp) strutfilmparameter
parameter determining strut-film boundary
Definition: globals.f90:18
real(dp) hi
film thickness at center
Definition: globals.f90:18
real(dp) int_reltol
relative tolerance of integrator
Definition: globals.f90:18
character(len=1024) viscositymodel
viscosity model (constant, fromFile)
Definition: globals.f90:10
real(dp) mu
viscosity
Definition: globals.f90:18
real(dp) cdp
disjoining pressure constant
Definition: globals.f90:18
real(dp) gr
growth rate
Definition: globals.f90:18
integer int_method
type of integration method (see MF in ODEPACK)
Definition: globals.f90:13
real(dp) bdp
disjoining pressure constant, set to zero for no DP
Definition: globals.f90:18
real(dp) ndp
disjoining pressure constant
Definition: globals.f90:18
real(dp) dstr
where strut starts in initial domain
Definition: globals.f90:18
real(dp) ae_tol
tolerance of algebraic equation solver
Definition: globals.f90:18
real(dp) timestep
timestep
Definition: globals.f90:49
integer its
number of outer integration outer time steps
Definition: globals.f90:29
real(dp) int_abstol
absolute tolerance of integrator
Definition: globals.f90:18
namespace with global variables
Definition: globals.f90:8
real(dp) rd
computational domain size (radius)
Definition: globals.f90:18
real(dp) initialtime
starting time
Definition: globals.f90:18
real(dp) hdp
disjoining pressure constant
Definition: globals.f90:18
integer maxts
maximum inner time steps between t and t+h
Definition: globals.f90:29
real(dp) mdp
disjoining pressure constant
Definition: globals.f90:18
real(dp) dr
mesh points spacing
Definition: globals.f90:18