~reducedmodelling/fluidity/ROM_Non-intrusive-ann

1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
!    Copyright (C) 2006 Imperial College London and others.
!    
!    Please see the AUTHORS file in the main source directory for a full list
!    of copyright holders.
!
!    Prof. C Pain
!    Applied Modelling and Computation Group
!    Department of Earth Science and Engineering
!    Imperial College London
!
!    amcgsoftware@imperial.ac.uk
!    
!    This library is free software; you can redistribute it and/or
!    modify it under the terms of the GNU Lesser General Public
!    License as published by the Free Software Foundation,
!    version 2.1 of the License.
!
!    This library is distributed in the hope that it will be useful,
!    but WITHOUT ANY WARRANTY; without even the implied warranty of
!    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
!    Lesser General Public License for more details.
!
!    You should have received a copy of the GNU Lesser General Public
!    License along with this library; if not, write to the Free Software
!    Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307
!    USA
#include "fdebug.h"

module multimaterial_module
  !! This module contains the options and material properties used
  !! when running FLUIDITY in the SOLIDITY mode
  use fldebug
  use state_module
  use fields
  use spud
  use fefields, only: compute_cv_mass
  use global_parameters, only: OPTION_PATH_LEN
  use field_priority_lists
  use field_options
  use diagnostic_fields_matrices
  use cv_upwind_values
  use equation_of_state, only: compressible_material_eos
  implicit none

  interface calculate_bulk_property
    module procedure calculate_bulk_scalar_property, calculate_bulk_vector_property, calculate_bulk_tensor_property
  end interface

  interface add_scaled_material_property
    module procedure add_scaled_material_property_scalar, add_scaled_material_property_vector, add_scaled_material_property_tensor
  end interface

  private
  public :: initialise_diagnostic_material_properties, &
            calculate_material_mass, calculate_bulk_material_pressure, &
            calculate_sum_material_volume_fractions, calculate_material_volume, &
            calculate_bulk_property, add_scaled_material_property, calculate_surfacetension, &
            calculate_diagnostic_material_volume_fraction, order_states_priority

contains

  subroutine calculate_surfacetension(state, surfacetension)
    ! calculates the surface tension in tensor form
    type(state_type), dimension(:), intent(inout) :: state
    type(tensor_field), intent(inout) :: surfacetension
    
    type(scalar_field), pointer :: volumefraction
    integer :: i, node, dimi, dimj, stat
    logical :: prognostic
    type(scalar_field) :: grad_mag, grad_mag2
    type(vector_field) :: gradient
    
    type(vector_field) :: normals
    logical, dimension(:), allocatable :: on_boundary
    type(vector_field), pointer :: x
    integer, dimension(2) :: shape_option
    integer, dimension(:), allocatable :: surface_ids
    
    real :: coeff, eq_angle
    real, dimension(surfacetension%dim(1), surfacetension%dim(2)) :: tensor
    
    if(size(state)==1) then
      FLExit("Don't know how to calculate a surface tension with only one material_phase.")
    end if
    
    x => extract_vector_field(state(1), "Coordinate")
    
    call allocate(gradient, surfacetension%dim(1), surfacetension%mesh, "Gradient")
    gradient%option_path = surfacetension%option_path
    
    call zero(surfacetension)
    
    do i = 1, size(state)
      volumefraction => extract_scalar_field(state(i), "MaterialVolumeFraction")
      
      prognostic = have_option(trim(volumefraction%option_path)//"/prognostic")
      
      if(prognostic.and.(.not.aliased(volumefraction))) then
        call get_option(trim(volumefraction%option_path)//&
                        "/prognostic/surface_tension/surface_tension_coefficient", &
                        coeff, default=0.0)
      
        call zero(gradient)
        
        call calculate_div_t_cv(state(i), gradient)
        
        ! the magnitude of the field gradient is a regularisation of the delta function
        ! indicating where the interface is
        grad_mag = magnitude(gradient)
        
        ! normalise the gradient
        do node = 1, node_count(surfacetension)
          if(node_val(grad_mag, node)>epsilon(0.0)) then
            call set(gradient, node, node_val(gradient, node)/node_val(grad_mag, node))
          else
            call set(gradient, node, spread(0.0, 1, gradient%dim))
            call set(grad_mag, node, 0.0)
          end if
        end do
        
        ! if we have an equilibrium contact angle then modify the gradient near the requested walls
        ! (note that grad_mag is unchanged)
        call get_option(trim(volumefraction%option_path)//&
                       "/prognostic/surface_tension/equilibrium_contact_angle", eq_angle, stat)
        if(stat==0) then
          call allocate(normals, mesh_dim(surfacetension), surfacetension%mesh, "NormalsToBoundary")
          call zero(normals)
          
          allocate(on_boundary(node_count(surfacetension)))
          on_boundary = .false.
        
          shape_option=option_shape(trim(volumefraction%option_path) // &
                 & "/prognostic/surface_tension/equilibrium_contact_angle/surface_ids")           
          allocate(surface_ids(1:shape_option(1)))
          call get_option(trim(volumefraction%option_path)//&
                 &"/prognostic/surface_tension/equilibrium_contact_angle/surface_ids", surface_ids)

          call calculate_boundary_normals(surfacetension%mesh, x, &
                                          normals, on_boundary, &
                                          surface_ids = surface_ids)
          
          do node = 1, node_count(surfacetension)
            if(on_boundary(node)) then
              call set(gradient, node, &
                      (node_val(normals, node)*cos(eq_angle)+node_val(gradient, node)*sin(eq_angle)))
            end if
          end do
          
          grad_mag2 = magnitude(gradient)
        
          ! renormalise the gradient
          do node = 1, node_count(surfacetension)
            if(node_val(grad_mag2, node)>epsilon(0.0)) then
              call set(gradient, node, node_val(gradient, node)/node_val(grad_mag2, node))
            else
              call set(gradient, node, spread(0.0, 1, gradient%dim))
            end if
          end do
          
          call deallocate(grad_mag2)
          
          deallocate(on_boundary)
          call deallocate(normals)
          deallocate(surface_ids)
          
        end if
                
        do node = 1, node_count(surfacetension)
          tensor = 0.0
          do dimi = 1, size(tensor,1)
            do dimj = 1, size(tensor,2)
              if(dimi==dimj) tensor(dimi,dimj) = coeff*node_val(grad_mag, node)
              tensor(dimi,dimj) = tensor(dimi,dimj) - &
                                  coeff*node_val(gradient, dimi, node)*node_val(gradient, dimj, node)*&
                                  node_val(grad_mag, node)
            end do
          end do
          
          call addto(surfacetension, node, tensor)
          
        end do
        
        call deallocate(grad_mag)
        
      end if
    
    end do
    
    call deallocate(gradient)
    
  end subroutine calculate_surfacetension

  subroutine initialise_diagnostic_material_properties(state)

    type(state_type), dimension(:), intent(inout) :: state

    !locals  
    integer :: stat, i
    type(scalar_field), pointer :: sfield
    logical :: prognostic

    do i = 1, size(state)

      sfield=>extract_scalar_field(state(i),'MaterialDensity',stat)
      if(stat==0) then
        prognostic=(have_option(trim(sfield%option_path)//'/prognostic'))
        if((.not.aliased(sfield)).and. prognostic) then
          call compressible_material_eos(state(i),materialdensity=sfield)
        end if
      end if

    end do

  end subroutine initialise_diagnostic_material_properties
  
  subroutine calculate_diagnostic_material_volume_fraction(state)

    type(state_type), dimension(:), intent(inout) :: state
    
    !locals
    type(scalar_field), pointer :: materialvolumefraction
    integer :: i, stat, diagnostic_count, diagnostic_state_index
    type(scalar_field) :: sumvolumefractions
    type(scalar_field), pointer :: sfield
    logical :: diagnostic

    ! How many diagnostic MaterialVolumeFraction fields do we have in state?
    ! Note that state contains all the submaterials of the current phase, including the phase itself.
    ! Therefore, if the only material is the phase itself, diagnostic_count should be 0. Otherwise,
    ! it should be 1.
    diagnostic_count = 0
    do i = 1, size(state)
       if(have_option(trim(state(i)%option_path)//"/scalar_field::MaterialVolumeFraction/diagnostic")) then
          diagnostic_count = diagnostic_count + 1
          ! Record the index of the state containing the diagnostic MaterialVolumeFraction field
          diagnostic_state_index = i 
       end if
    end do

    if(diagnostic_count>1) then
      ewrite(0,*) diagnostic_count, ' diagnostic MaterialVolumeFractions'
      FLExit("Only one diagnostic MaterialVolumeFraction permitted.")
    end if

    if(diagnostic_count==1) then
      ! Extract the diagnostic volume fraction
      materialvolumefraction => extract_scalar_field(state(diagnostic_state_index), 'MaterialVolumeFraction')
      
      call allocate(sumvolumefractions, materialvolumefraction%mesh, 'Sum of volume fractions')
      call zero(sumvolumefractions)
      
      do i = 1,size(state)
        sfield=>extract_scalar_field(state(i),'MaterialVolumeFraction',stat)
        diagnostic=(have_option(trim(sfield%option_path)//'/diagnostic'))
        if ( (stat==0).and.(.not. aliased(sfield)).and.(.not.diagnostic)) then
          call addto(sumvolumefractions, sfield)
        end if
      end do
      
      call set(materialvolumefraction, 1.0)
      call addto(materialvolumefraction, sumvolumefractions, -1.0)
      call deallocate(sumvolumefractions)
    end if

  end subroutine calculate_diagnostic_material_volume_fraction
  
  subroutine order_states_priority(state, state_order)
    type(state_type), dimension(:), intent(inout) :: state
    integer, dimension(:), intent(inout) :: state_order
    
    type(scalar_field), pointer :: volumefraction
    logical, dimension(size(state)) :: priority_states
    integer, dimension(size(state)) :: state_priorities

    integer :: i, p, f

    assert(size(state_order)==size(state))

    priority_states = .false.
    state_priorities = 0
    do i = 1, size(state)
      volumefraction=>extract_scalar_field(state(i), "MaterialVolumeFraction")

      if(have_option(trim(volumefraction%option_path)//"/prognostic/priority")) then
        call get_option(trim(volumefraction%option_path)//"/prognostic/priority", state_priorities(i))
        priority_states(i) = .true.
      end if
    end do

    do i = 1, size(state)
      if(.not.priority_states(i)) then
        state_priorities(i) = minval(state_priorities)-1
      end if
    end do

    ! now work out the right order
    f = 0
    state_order = 0
    do p = maxval(state_priorities), minval(state_priorities), -1
      do i=1, size(state)
        if(state_priorities(i)==p) then
          f = f + 1
          state_order(f) = i
        end if
      end do
    end do
    assert(all(state_order>0))

  end subroutine order_states_priority

  subroutine calculate_bulk_scalar_property(state,bulkfield,materialname,momentum_diagnostic)

    type(state_type), dimension(:), intent(inout) :: state
    type(scalar_field), intent(inout) :: bulkfield
    character(len=*), intent(in) :: materialname
    logical, intent(in), optional ::momentum_diagnostic

    !locals
    integer :: i, stat
    type(scalar_field), pointer :: sfield
 
    integer, dimension(size(state)) :: state_order
    type(scalar_field) :: sumvolumefractionsbound

    ewrite(1,*) 'In calculate_bulk_scalar_property'

    call zero(bulkfield)

    call order_states_priority(state, state_order)

    call allocate(sumvolumefractionsbound, bulkfield%mesh, "SumMaterialVolumeFractionsBound")
    call set(sumvolumefractionsbound, 1.0)
    
    do i = 1, size(state)
  
      ewrite(2,*) 'Considering state: ', state(state_order(i))%name
      sfield => extract_scalar_field(state(state_order(i)), trim(materialname), stat)
      if(stat==0) then
        call add_scaled_material_property(state(state_order(i)), bulkfield, sfield, &
                                          sumvolumefractionsbound=sumvolumefractionsbound, &
                                          momentum_diagnostic=momentum_diagnostic)
      end if

    end do
    
    call deallocate(sumvolumefractionsbound)
    
  end subroutine calculate_bulk_scalar_property

  subroutine calculate_bulk_vector_property(state,bulkfield,materialname,momentum_diagnostic)

    type(state_type), dimension(:), intent(inout) :: state
    type(vector_field), intent(inout) :: bulkfield
    character(len=*), intent(in) :: materialname
    logical, intent(in), optional :: momentum_diagnostic

    !locals
    integer :: i, stat
    type(vector_field), pointer :: vfield
 
    integer, dimension(size(state)) :: state_order
    type(scalar_field) :: sumvolumefractionsbound

    ewrite(1,*) 'In calculate_bulk_vector_property'

    call zero(bulkfield)
    
    call order_states_priority(state, state_order)
    
    call allocate(sumvolumefractionsbound, bulkfield%mesh, "SumMaterialVolumeFractionsBound")
    call set(sumvolumefractionsbound, 1.0)

    do i = 1, size(state)
  
      ewrite(2,*) 'Considering state: ', state(state_order(i))%name
      vfield => extract_vector_field(state(state_order(i)), trim(materialname), stat)
      if(stat==0) then
        call add_scaled_material_property(state(state_order(i)), bulkfield, vfield, &
                                          sumvolumefractionsbound=sumvolumefractionsbound, &
                                          momentum_diagnostic=momentum_diagnostic)
      end if

    end do
    
    call deallocate(sumvolumefractionsbound)

  end subroutine calculate_bulk_vector_property

  subroutine calculate_bulk_tensor_property(state,bulkfield,materialname,momentum_diagnostic)

    type(state_type), dimension(:), intent(inout) :: state
    type(tensor_field), intent(inout) :: bulkfield
    character(len=*), intent(in) :: materialname
    logical, intent(in), optional :: momentum_diagnostic

    !locals
    integer :: i, stat
    type(tensor_field), pointer :: tfield
 
    integer, dimension(size(state)) :: state_order
    type(scalar_field) :: sumvolumefractionsbound

    ewrite(1,*) 'In calculate_bulk_tensor_property'

    call zero(bulkfield)
    
    call order_states_priority(state, state_order)
    
    call allocate(sumvolumefractionsbound, bulkfield%mesh, "SumMaterialVolumeFractionsBound")
    call set(sumvolumefractionsbound, 1.0)

    do i = 1, size(state)
  
      ewrite(2,*) 'Considering state: ', state(state_order(i))%name
      tfield => extract_tensor_field(state(state_order(i)), trim(materialname), stat)
      if(stat==0) then
        call add_scaled_material_property(state(state_order(i)), bulkfield, tfield, &
                                          sumvolumefractionsbound=sumvolumefractionsbound, &
                                          momentum_diagnostic=momentum_diagnostic)
      end if
    
    end do

    call deallocate(sumvolumefractionsbound)
    
  end subroutine calculate_bulk_tensor_property

  subroutine get_scalable_volume_fraction(scaledvfrac, state, sumvolumefractionsbound, momentum_diagnostic)

    type(scalar_field), intent(inout) :: scaledvfrac
    type(state_type), intent(inout) :: state
    type(scalar_field), intent(inout), optional :: sumvolumefractionsbound
    logical, intent(in), optional :: momentum_diagnostic

    type(scalar_field), pointer :: volumefraction, oldvolumefraction
    type(vector_field), pointer :: velocity
    type(scalar_field) :: remapvfrac

    integer :: stat
    real :: theta
    
    logical :: cap
    real:: u_cap_val, l_cap_val

    volumefraction => extract_scalar_field(state, 'MaterialVolumeFraction')
    
    call remap_field(volumefraction, scaledvfrac)
          
    if(present_and_true(momentum_diagnostic)) then
      velocity => extract_vector_field(state, 'Velocity', stat=stat)
      if(stat==0) then
        call get_option(trim(velocity%option_path)//'/prognostic/temporal_discretisation/theta', &
                        theta, stat)
        if(stat==0) then
          call allocate(remapvfrac, scaledvfrac%mesh, "RemppedMaterialVolumeFraction")
          
          oldvolumefraction => extract_scalar_field(state, 'OldMaterialVolumeFraction')
          call remap_field(oldvolumefraction, remapvfrac)
          
          call scale(scaledvfrac, theta)
          call addto(scaledvfrac, remapvfrac, (1.-theta))
          
          call deallocate(remapvfrac)
        end if
      end if
    end if
          
    cap = (have_option(trim(complete_field_path(volumefraction%option_path))//"/cap_values"))
          
    if(cap) then
      ! this capping takes care of under or overshoots in this volume fraction individually
      ! these will have typically occurred during advection
            
      call get_option(trim(complete_field_path(volumefraction%option_path))//"/cap_values/upper_cap", &
                      u_cap_val, default=huge(0.0)*epsilon(0.0))
      call get_option(trim(complete_field_path(volumefraction%option_path))//"/cap_values/lower_cap", &
                      l_cap_val, default=-huge(0.0)*epsilon(0.0))
            
      call bound(scaledvfrac, l_cap_val, u_cap_val)
          
      if(present(sumvolumefractionsbound)) then
        assert(sumvolumefractionsbound%mesh==scaledvfrac%mesh)
        ! this capping takes care of overlapping volume fractions
        call bound(scaledvfrac, upper_bound=sumvolumefractionsbound) 
        call addto(sumvolumefractionsbound, scaledvfrac, scale=-1.0)
        ewrite_minmax(sumvolumefractionsbound)
      end if
          
    end if
    ewrite_minmax(scaledvfrac)
 
  end subroutine get_scalable_volume_fraction

  subroutine add_scaled_material_property_scalar(state,bulkfield,field,sumvolumefractionsbound,momentum_diagnostic)

    type(state_type), intent(inout) :: state
    type(scalar_field), intent(inout) :: bulkfield, field
    type(scalar_field), intent(inout), optional :: sumvolumefractionsbound
    logical, intent(in), optional :: momentum_diagnostic

    !locals
    type(scalar_field) :: scaledvfrac
    type(scalar_field) :: tempfield
    
    call allocate(tempfield, bulkfield%mesh, "Temp"//trim(bulkfield%name))
    call allocate(scaledvfrac, bulkfield%mesh, "ScaledMaterialVolumeFraction")

    call get_scalable_volume_fraction(scaledvfrac, state, &
                                      sumvolumefractionsbound=sumvolumefractionsbound, &
                                      momentum_diagnostic=momentum_diagnostic)

    call remap_field(field, tempfield)
    call scale(tempfield, scaledvfrac)
    call addto(bulkfield, tempfield)

    call deallocate(tempfield)
    call deallocate(scaledvfrac)

  end subroutine add_scaled_material_property_scalar

  subroutine add_scaled_material_property_vector(state,bulkfield,field,sumvolumefractionsbound,momentum_diagnostic)

    type(state_type), intent(inout) :: state
    type(vector_field), intent(inout) :: bulkfield, field
    type(scalar_field), intent(inout), optional :: sumvolumefractionsbound
    logical, intent(in), optional :: momentum_diagnostic

    !locals
    type(scalar_field) :: scaledvfrac
    type(vector_field) :: tempfield

    call allocate(tempfield, bulkfield%dim, bulkfield%mesh, "Temp"//trim(bulkfield%name))
    call allocate(scaledvfrac, bulkfield%mesh, "ScaledMaterialVolumeFraction")
    
    call get_scalable_volume_fraction(scaledvfrac, state, &
                                      sumvolumefractionsbound=sumvolumefractionsbound, &
                                      momentum_diagnostic=momentum_diagnostic)
          
    call remap_field(field, tempfield)
    call scale(tempfield, scaledvfrac)
    call addto(bulkfield, tempfield)
          
    call deallocate(tempfield)
    call deallocate(scaledvfrac)

  end subroutine add_scaled_material_property_vector

  subroutine add_scaled_material_property_tensor(state,bulkfield,field,sumvolumefractionsbound,momentum_diagnostic)

    type(state_type), intent(inout) :: state
    type(tensor_field), intent(inout) :: bulkfield, field
    type(scalar_field), intent(inout), optional :: sumvolumefractionsbound
    logical, intent(in), optional :: momentum_diagnostic

    !locals
    type(scalar_field) :: scaledvfrac
    type(tensor_field) :: tempfield

    call allocate(tempfield, bulkfield%mesh, "Temp"//trim(bulkfield%name))
    call allocate(scaledvfrac, bulkfield%mesh, "ScaledMaterialVolumeFraction")

    call get_scalable_volume_fraction(scaledvfrac, state, &
                                      sumvolumefractionsbound=sumvolumefractionsbound, &
                                      momentum_diagnostic=momentum_diagnostic)
    
    call remap_field(field, tempfield)
    call scale(tempfield, scaledvfrac)
    call addto(bulkfield, tempfield)
          
    call deallocate(tempfield)
    call deallocate(scaledvfrac)

  end subroutine add_scaled_material_property_tensor

  subroutine calculate_material_volume(state, materialvolume)

    type(state_type), intent(in) :: state
    type(scalar_field), intent(inout) :: materialvolume
  
    ! local
    type(scalar_field) :: cvmass
    type(scalar_field), pointer :: volumefraction
    type(vector_field), pointer :: coordinates
  
    coordinates=>extract_vector_field(state, "Coordinate")
  
    call allocate(cvmass, materialvolume%mesh, "CV mass")
    call zero(cvmass)
  
    call compute_cv_mass(coordinates, cvmass)
  
    volumefraction=>extract_scalar_field(state,"MaterialVolumeFraction")
  
    materialvolume%val=volumefraction%val*cvmass%val
  
    call deallocate(cvmass)

  end subroutine calculate_material_volume

  subroutine calculate_material_mass(state, materialmass)

    type(state_type), intent(in) :: state
    type(scalar_field), intent(inout) :: materialmass
  
    ! local
    integer :: stat
    type(scalar_field) :: cvmass
    type(scalar_field), pointer :: volumefraction, materialdensity
    type(vector_field), pointer :: coordinates
    real :: rho_0
  
    coordinates=>extract_vector_field(state, "Coordinate")
  
    call allocate(cvmass, materialmass%mesh, "CV mass")
    call zero(cvmass)
  
    call compute_cv_mass(coordinates, cvmass)

    volumefraction=>extract_scalar_field(state,"MaterialVolumeFraction")
    call set(materialmass, volumefraction)
    call scale(materialmass, cvmass)
  
    materialdensity=>extract_scalar_field(state,"MaterialDensity", stat=stat)
    if(stat==0) then
      call scale(materialmass, materialdensity)
    else
      call get_option("/material_phase::"//trim(state%name)&
                      //"/equation_of_state/fluids/linear/reference_density", rho_0)
      call scale(materialmass, rho_0)
    end if
  
    call deallocate(cvmass)

  end subroutine calculate_material_mass

  subroutine calculate_bulk_material_pressure(state,bulkmaterialpressure)

    type(state_type), dimension(:), intent(inout) :: state
    type(scalar_field) :: bulkmaterialpressure

    !locals
    integer :: i, stat
    type(scalar_field), pointer :: volumefraction
    type(scalar_field) :: materialpressure

    call zero(bulkmaterialpressure)

    call allocate(materialpressure, bulkmaterialpressure%mesh, "TempBulkPressure")

    do i = 1, size(state)

      volumefraction => extract_scalar_field(state(i), 'MaterialVolumeFraction', stat)

      if (stat==0) then

         call compressible_material_eos(state(i), materialpressure=materialpressure)

         bulkmaterialpressure%val=bulkmaterialpressure%val+volumefraction%val*materialpressure%val

      end if

    end do

    call deallocate(materialpressure)

  end subroutine calculate_bulk_material_pressure
  
  subroutine calculate_sum_material_volume_fractions(state,sumvolumefractions)

    type(state_type), dimension(:), intent(inout) :: state
    type(scalar_field), intent(inout) :: sumvolumefractions

    !locals
    integer :: i, stat
    type(scalar_field), pointer :: sfield
    logical :: prognostic, diagnostic, prescribed

    diagnostic = have_option(trim(sumvolumefractions%option_path)//"/diagnostic")
    if(.not.diagnostic) return

    call zero(sumvolumefractions)

    do i = 1,size(state)
      sfield=>extract_scalar_field(state(i),'MaterialVolumeFraction',stat)
      if(stat==0) then
        prognostic = have_option(trim(sfield%option_path)//"/prognostic")
        prescribed = have_option(trim(sfield%option_path)//"/prescribed")
        if ((.not.aliased(sfield)).and.(prognostic.or.prescribed)) then
          call addto(sumvolumefractions, sfield)
        end if
      end if
    end do

  end subroutine calculate_sum_material_volume_fractions

end module multimaterial_module