MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
InOut.f90
Go to the documentation of this file.
1 
8 module inout
9  implicit none
10  private
11  public input,output,print_header
12 contains
13 !********************************BEGINNING*************************************
17 subroutine input()
18  use constants
19  use globals
20  use physicalproperties
21  use ioutils
22  use fson
23  use fson_value_m, only: fson_value_get
24  type(fson_value), pointer :: json_data
25  character(len=1024) :: strval
26  character(len=99) :: after_foaming,after_foaming0='after_foaming.txt'
27  character(len=99) :: bg_res='../../foamExpansion/results/bubbleGrowth/'
28  character(len=99) :: qmom0D_res='../../foamExpansion/results/CFD0D/'
29  character(len=99) :: qmom3D_res='../../foamExpansion/results/CFD3D/'
30  real(dp) :: matr(7)
31  integer :: fi
32  ! Read input parameters
33  json_data => fson_parse("../inputs/foamAging.json")
34  call fson_get(json_data, "numerics.timeStart", tbeg)
35  call fson_get(json_data, "numerics.timeEnd", tend)
36  call fson_get(json_data, "numerics.numberOfOutputs", nroutputs)
37  call fson_get(json_data, "numerics.wallDiscretization", divwall)
38  call fson_get(json_data, "numerics.cellDiscretization", divcell)
39  call fson_get(json_data, "numerics.sheetDiscretization", divsheet)
40  call fson_get(json_data, "foamCondition.foamHalfThickness", dfoam)
41  call fson_get(json_data, "foamCondition.inProtectiveSheet", sheet)
42  if (sheet) then
43  call fson_get(&
44  json_data, "foamCondition.sheetThickness", dsheet)
45  endif
46  call fson_get(json_data, "foamCondition.agingTemperature", temp)
47  call fson_get(json_data, "foamCondition.conductivityTemperature", temp_cond)
48  call fson_get(json_data, "foamCondition.initialPressure", pressure)
49  call fson_get(json_data, "foamCondition.boundaryPressure.O2", pbg(1))
50  call fson_get(json_data, "foamCondition.boundaryPressure.N2", pbg(2))
51  call fson_get(json_data, "foamCondition.boundaryPressure.CO2", pbg(3))
52  call fson_get(&
53  json_data, "foamCondition.boundaryPressure.Cyclopentane", pbg(4))
54  call fson_get(json_data, "sourceOfProperty.gasComposition", strval)
55  if (strval=="BubbleGrowth") then
56  after_foaming=trim(adjustl(bg_res))//trim(adjustl(after_foaming0))
57  open(unit=newunit(fi),file=after_foaming)
58  read(fi,*)
59  read(fi,*) matr(1:4)
60  xg(1)=0 ! no oxygen
61  xg(2)=0 ! no nitrogen
62  xg(3)=matr(4)
63  xg(4)=matr(3)
64  xg=xg/sum(xg)
65  close(fi)
66  elseif (strval=="Qmom0D") then
67  after_foaming=trim(adjustl(qmom0d_res))//trim(adjustl(after_foaming0))
68  open(unit=newunit(fi),file=after_foaming)
69  read(fi,*)
70  read(fi,*) matr(1:4)
71  xg(1)=0 ! no oxygen
72  xg(2)=0 ! no nitrogen
73  xg(3)=matr(4)
74  xg(4)=matr(3)
75  xg=xg/sum(xg)
76  close(fi)
77  elseif (strval=="Qmom3D") then
78  after_foaming=trim(adjustl(qmom3d_res))//trim(adjustl(after_foaming0))
79  open(unit=newunit(fi),file=after_foaming)
80  read(fi,*)
81  read(fi,*) matr(1:4)
82  xg(1)=0 ! no oxygen
83  xg(2)=0 ! no nitrogen
84  xg(3)=matr(4)
85  xg(4)=matr(3)
86  xg=xg/sum(xg)
87  close(fi)
88  elseif (strval=="DirectInput") then
89  call fson_get(json_data, "foamCondition.initialComposition.O2", xg(1))
90  call fson_get(json_data, "foamCondition.initialComposition.N2", xg(2))
91  call fson_get(json_data, "foamCondition.initialComposition.CO2", xg(3))
92  call fson_get(&
93  json_data, "foamCondition.initialComposition.Cyclopentane", xg(4))
94  xg=xg/sum(xg)
95  else
96  write(*,*) 'unknown source for gas composition'
97  stop
98  endif
99  call fson_get(json_data, "sourceOfProperty.foamDensity", strval)
100  if (strval=="BubbleGrowth") then
101  after_foaming=trim(adjustl(bg_res))//trim(adjustl(after_foaming0))
102  open(unit=newunit(fi),file=after_foaming)
103  read(fi,*)
104  read(fi,*) matr(1:4)
105  rhof=matr(1)
106  close(fi)
107  elseif (strval=="Qmom0D") then
108  after_foaming=trim(adjustl(qmom0d_res))//trim(adjustl(after_foaming0))
109  open(unit=newunit(fi),file=after_foaming)
110  read(fi,*)
111  read(fi,*) matr(1:4)
112  rhof=matr(1)
113  close(fi)
114  elseif (strval=="Qmom3D") then
115  after_foaming=trim(adjustl(qmom3d_res))//trim(adjustl(after_foaming0))
116  open(unit=newunit(fi),file=after_foaming)
117  read(fi,*)
118  read(fi,*) matr(1:4)
119  rhof=matr(1)
120  close(fi)
121  elseif (strval=="DirectInput") then
122  call fson_get(json_data, "morphology.foamDensity", rhof)
123  else
124  write(*,*) 'unknown source for foam density'
125  stop
126  endif
127  call fson_get(json_data, "sourceOfProperty.cellSize", strval)
128  if (strval=="BubbleGrowth") then
129  after_foaming=trim(adjustl(bg_res))//trim(adjustl(after_foaming0))
130  open(unit=newunit(fi),file=after_foaming)
131  read(fi,*)
132  read(fi,*) matr(1:4)
133  dcell=matr(2)
134  close(fi)
135  elseif (strval=="Qmom0D") then
136  after_foaming=trim(adjustl(qmom0d_res))//trim(adjustl(after_foaming0))
137  open(unit=newunit(fi),file=after_foaming)
138  read(fi,*)
139  read(fi,*) matr(1:4)
140  dcell=matr(2)
141  close(fi)
142  elseif (strval=="Qmom3D") then
143  after_foaming=trim(adjustl(qmom3d_res))//trim(adjustl(after_foaming0))
144  open(unit=newunit(fi),file=after_foaming)
145  read(fi,*)
146  read(fi,*) matr(1:4)
147  dcell=matr(2)
148  close(fi)
149  elseif (strval=="DirectInput") then
150  call fson_get(json_data, "morphology.cellSize", dcell)
151  else
152  write(*,*) 'unknown source for cell size'
153  stop
154  endif
155  call fson_get(json_data, "sourceOfProperty.strutContent", strval)
156  if (strval=="StrutContent") then
157  call strutcontent(fstrut,rhof)
158  elseif (strval=="DirectInput") then
159  call fson_get(json_data, "morphology.strutContent", fstrut)
160  else
161  write(*,*) 'unknown source for strut content'
162  stop
163  endif
164  call fson_get(json_data, "sourceOfProperty.wallThickness", strval)
165  if (strval=="DirectInput") then
166  call fson_get(json_data, "morphology.wallThickness", dwall)
167  else
168  write(*,*) 'unknown source for wall thickness'
169  stop
170  endif
171  call fson_get(json_data, "physicalProperties.polymerDensity", rhop)
172  call fson_get(json_data, &
173  "physicalProperties.foam.solubilityModel.O2", strval)
174  if (strval=="constant") then
175  solmodel(1)=0
176  elseif ( strval=="modena" ) then
177  solmodel(1)=1
178  else
179  print*, "Solubility model must be constant or modena"
180  endif
181  call fson_get(json_data, &
182  "physicalProperties.foam.solubilityModel.N2", strval)
183  if (strval=="constant") then
184  solmodel(2)=0
185  elseif ( strval=="modena" ) then
186  solmodel(2)=1
187  else
188  print*, "Solubility model must be constant or modena"
189  endif
190  call fson_get(json_data, &
191  "physicalProperties.foam.solubilityModel.CO2", strval)
192  if (strval=="constant") then
193  solmodel(3)=0
194  elseif ( strval=="modena" ) then
195  solmodel(3)=1
196  else
197  print*, "Solubility model must be constant or modena"
198  endif
199  call fson_get(json_data, &
200  "physicalProperties.foam.solubilityModel.Cyclopentane", strval)
201  if (strval=="constant") then
202  solmodel(4)=0
203  elseif ( strval=="modena" ) then
204  solmodel(4)=1
205  else
206  print*, "Solubility model must be constant or modena"
207  endif
208  if (solmodel(1)==0) then
209  call fson_get(json_data, "physicalProperties.foam.solubility.O2", sg(1))
210  endif
211  if (solmodel(2)==0) then
212  call fson_get(json_data, "physicalProperties.foam.solubility.N2", sg(2))
213  endif
214  if (solmodel(3)==0) then
215  call fson_get(&
216  json_data, "physicalProperties.foam.solubility.CO2", sg(3))
217  endif
218  if (solmodel(4)==0) then
219  call fson_get(&
220  json_data, "physicalProperties.foam.solubility.Cyclopentane", sg(4))
221  endif
222  call fson_get(json_data, &
223  "physicalProperties.foam.diffusivityModel.O2", strval)
224  if (strval=="constant") then
225  diffmodel(1)=0
226  elseif ( strval=="modena" ) then
227  diffmodel(1)=1
228  else
229  print*, "Diffusivity model must be constant or modena"
230  endif
231  call fson_get(json_data, &
232  "physicalProperties.foam.diffusivityModel.N2", strval)
233  if (strval=="constant") then
234  diffmodel(2)=0
235  elseif ( strval=="modena" ) then
236  diffmodel(2)=1
237  else
238  print*, "Diffusivity model must be constant or modena"
239  endif
240  call fson_get(json_data, &
241  "physicalProperties.foam.diffusivityModel.CO2",strval)
242  if (strval=="constant") then
243  diffmodel(3)=0
244  elseif ( strval=="modena" ) then
245  diffmodel(3)=1
246  else
247  print*, "Diffusivity model must be constant or modena"
248  endif
249  call fson_get(json_data, &
250  "physicalProperties.foam.diffusivityModel.Cyclopentane", strval)
251  if (strval=="constant") then
252  diffmodel(4)=0
253  elseif ( strval=="modena" ) then
254  diffmodel(4)=1
255  else
256  print*, "Diffusivity model must be constant or modena"
257  endif
258  if (diffmodel(1)==0) then
259  call fson_get(json_data, &
260  "physicalProperties.foam.diffusivity.O2", dg(1))
261  endif
262  if (diffmodel(2)==0) then
263  call fson_get(json_data, &
264  "physicalProperties.foam.diffusivity.N2", dg(2))
265  endif
266  if (diffmodel(3)==0) then
267  call fson_get(json_data, &
268  "physicalProperties.foam.diffusivity.CO2", dg(3))
269  endif
270  if (diffmodel(4)==0) then
271  call fson_get(json_data, &
272  "physicalProperties.foam.diffusivity.Cyclopentane", dg(4))
273  endif
274  if (sheet) then
275  call fson_get(json_data, &
276  "physicalProperties.sheet.solubility.O2", sheetsg(1))
277  call fson_get(json_data, &
278  "physicalProperties.sheet.solubility.N2", sheetsg(2))
279  call fson_get(json_data, &
280  "physicalProperties.sheet.solubility.CO2", sheetsg(3))
281  call fson_get(json_data, &
282  "physicalProperties.sheet.solubility.Cyclopentane", sheetsg(4))
283  call fson_get(json_data, &
284  "physicalProperties.sheet.diffusivity.O2", sheetdg(1))
285  call fson_get(json_data, &
286  "physicalProperties.sheet.diffusivity.N2", sheetdg(2))
287  call fson_get(json_data, &
288  "physicalProperties.sheet.diffusivity.CO2", sheetdg(3))
289  call fson_get(json_data, &
290  "physicalProperties.sheet.diffusivity.Cyclopentane", sheetdg(4))
291  endif
292  call fson_destroy(json_data)
293 end subroutine input
294 !***********************************END****************************************
295 
296 
297 !********************************BEGINNING*************************************
302 subroutine output(iprof, time, ystate, neq, pp)
303  use constants
304  use globals
305  use model, only: ngas,dz,mor,nfv,sol
306  integer, intent(in) :: iprof
307  integer, intent(in) :: neq
308  real(dp), intent(in) :: time
309  real(dp), intent(in) :: ystate(:)
310  real(dp), intent(out) :: pp(:)
311  integer :: i, j
312  integer :: spp
313  real(dp) :: pos
314 
315  character(len=1) :: name_1 ! one character
316  character(len=2) :: name_2 ! two characters
317  character(len=3) :: name_3 ! three characters
318  character(len=4) :: name_f ! final name of file
319 
320  if (iprof < 10) then
321  write(name_1,'(I1)') iprof
322  name_f = '000' // name_1
323  elseif (iprof >= 10 .and. iprof < 100) then
324  write(name_2,'(I2)') iprof
325  name_f = '00' // name_2
326  elseif (iprof >= 100 .and. iprof < 1000) then
327  write(name_3,'(I3)') iprof
328  name_f = '0' // name_3
329  else
330  write(name_f,'(I4)') iprof
331  endif
332  open(unit=11,file='H2perm_'//trim(name_f)//'.dat')
333  open(unit=12,file='ppar_'//trim(name_f)//'.dat')
334 
335  ! profiles
336  do i = 1, neq/ngas
337  write (11,100) time/(3600*24),pos,&
338  ystate(ngas*(i-1)+1)*sol(ngas*(i-1)+1),&
339  ystate(ngas*(i-1)+2)*sol(ngas*(i-1)+2),&
340  ystate(ngas*(i-1)+3)*sol(ngas*(i-1)+3),&
341  ystate(ngas*(i-1)+4)*sol(ngas*(i-1)+4)
342  enddo
343  pp=0
344  spp=0
345  do i = 1,nfv
346  if (i==1) then
347  pos=dz(1)
348  else
349  pos=pos+(dz(i-1)+dz(i))/2
350  endif
351  if (mor(i)==1) then
352  do j=1,ngas
353  pp(j)=pp(j)+ystate(ngas*(i-1)+j)*rg*temp
354  enddo
355  spp=spp+1
356  write (12,101) time/(3600*24),pos,ystate(ngas*(i-1)+1)*rg*temp,&
357  ystate(ngas*(i-1)+2)*rg*temp,ystate(ngas*(i-1)+3)*rg*temp,&
358  ystate(ngas*(i-1)+4)*rg*temp
359  endif
360  enddo
361  pp=pp/spp
362  close(11)
363  close(12)
364 100 format (f8.2,2x,es23.8e3,2x,es23.8e3,2x,es23.8e3,2x,es23.8e3,2x,es23.8e3)
365 101 format (f8.2,2x,es23.8e3,2x,es23.8e3,2x,es23.8e3,2x,es23.8e3,2x,es23.8e3)
366 end subroutine output
367 !***********************************END****************************************
368 
369 
370 !********************************BEGINNING*************************************
374 subroutine print_header
376  use model, only: nfv
377  print*, 'Foam:'
378  write(*,'(A30,EN12.3,1x,A)') 'foam density:',rhof,'kg/m3'
379  write(*,'(A30,EN12.3)') 'porosity:',eps
380  write(*,'(A30,EN12.3)') 'strut content:',fstrut
381  write(*,'(A30,EN12.3,1x,A)') 'cell size:',dcell,'m'
382  write(*,'(A30,EN12.3,1x,A)') 'wall thickness:',dwall,'m'
383  if (sheet) then
384  write(*,'(A30,EN12.3,1x,A)') 'sheet thickness:',dsheet,'m'
385  endif
386  write(*,'(A30,EN12.3,1x,A)') 'foam thickness:',dfoam,'m'
387  write(*,'(A30,I12)') 'number of cells:',ncell
388  write(*,'(A30,EN12.3,1x,A)') 'aging temperature:',temp,'K'
389  write(*,'(A30,EN12.3,1x,A)') 'conductivity temperature:',temp_cond,'K'
390  print*, 'Physical properties:'
391  write(*,'(A30,EN12.3,1x,A)') 'polymer density:',rhop,'kg/m3'
392  write(*,'(A30,EN12.3,1x,A)') 'diffusivity in gas:',dgas,'m2/s'
393  write(*,'(A30,EN12.3,1x,A)') 'O2 diffusivity:',dg(1),'m2/s'
394  write(*,'(A30,EN12.3,1x,A)') 'N2 diffusivity:',dg(2),'m2/s'
395  write(*,'(A30,EN12.3,1x,A)') 'CO2 diffusivity:',dg(3),'m2/s'
396  write(*,'(A30,EN12.3,1x,A)') 'pentane diffusivity:',dg(4),'m2/s'
397  write(*,'(A30,EN12.3,1x,A)') 'O2 solubility:',sg(1),'g/g/bar'
398  write(*,'(A30,EN12.3,1x,A)') 'N2 solubility:',sg(2),'g/g/bar'
399  write(*,'(A30,EN12.3,1x,A)') 'CO2 solubility:',sg(3),'g/g/bar'
400  write(*,'(A30,EN12.3,1x,A)') 'pentane solubility:',sg(4),'g/g/bar'
401  write(*,'(A30,EN12.3,1x,A)') 'O2 permeability:',sg(1)*dg(1),'m2/s*g/g/bar'
402  write(*,'(A30,EN12.3,1x,A)') 'N2 permeability:',sg(2)*dg(2),'m2/s*g/g/bar'
403  write(*,'(A30,EN12.3,1x,A)') 'CO2 permeability:',sg(3)*dg(3),'m2/s*g/g/bar'
404  write(*,'(A30,EN12.3,1x,A)') 'pentane permeability:',sg(4)*dg(4),&
405  'm2/s*g/g/bar'
406  print*, 'Numerics:'
407  write(*,'(A30,I12)') 'finite volumes in wall:',divwall
408  write(*,'(A30,I12)') 'finite volumes in cell:',divcell
409  if (sheet) then
410  write(*,'(A30,I12)') 'finite volumes in sheet:',divsheet
411  endif
412  write(*,'(A30,I12)') 'finite volumes in total:',nfv
413  write(*,'(A30,EN12.3,1x,A)') 'initial time:',tbeg,'s'
414  write(*,'(A30,EN12.3,1x,A)') 'end time:',tend,'s'
415 end subroutine print_header
416 !***********************************END****************************************
417 end module inout
real(dp) dcell
cell size
Definition: globals.f90:20
integer divwall
number of grid points in wall
Definition: globals.f90:12
integer, dimension(:), allocatable solmodel
solubility model
Definition: globals.f90:16
real(dp), dimension(:), allocatable xg
molar fractions of gases
Definition: globals.f90:36
integer ngas
number of dissolved gases
Definition: globals.f90:29
real(dp), dimension(:), allocatable sheetsg
gas solubility in sheet
Definition: globals.f90:34
integer divcell
number of grid points in cell
Definition: globals.f90:13
real(dp) rhop
polymer density
Definition: globals.f90:49
real(dp) temp
temperature (K)
Definition: globals.f90:49
real(dp) tbeg
time at the end
Definition: globals.f90:19
real(dp) eps
foam porosity
Definition: globals.f90:30
integer, dimension(:), allocatable diffmodel
diffusivity model
Definition: globals.f90:17
real(dp), dimension(:), allocatable pbg
boundary conditions
Definition: globals.f90:37
real(dp), dimension(:), allocatable dz
spatial discretization
Definition: globals.f90:106
real(dp), dimension(:), allocatable pressure
partial pressure(Pa)
Definition: globals.f90:106
real(dp), dimension(:), allocatable dg
gas diffusiviy in polymer
Definition: globals.f90:33
real(dp) dfoam
foam thickness
Definition: globals.f90:22
real(dp) rhof
foam density
Definition: globals.f90:28
integer nroutputs
number of outer endme steps
Definition: globals.f90:11
real(dp) fstrut
strut content
Definition: globals.f90:27
real(dp), dimension(:), allocatable sheetdg
gas diffusivity in sheet
Definition: globals.f90:35
integer divsheet
number of grid points in sheet
Definition: globals.f90:14
real(dp), dimension(:), allocatable sg
gas solubility in polymer
Definition: globals.f90:32
integer ncell
number of cells
Definition: globals.f90:15
real(dp) dsheet
sheet thickness
Definition: globals.f90:23
real(dp) temp_cond
temperature of conductivity measurements
Definition: globals.f90:25
namespace with global variables
Definition: globals.f90:8
logical sheet
determines if sheet is used on the outside of the foam
Definition: globals.f90:10
real(dp) dgas
gas diffusivity in gasphase
Definition: globals.f90:31
real(dp) dwall
wall thickness
Definition: globals.f90:21
double tend
end time, s