richards_common.F90       coverage:  100.00 %func     76.53 %block


     1) module Richards_Common_module
     2) 
     3)   use Richards_Aux_module
     4)   use Global_Aux_module
     5)   
     6)   use PFLOTRAN_Constants_module
     7) 
     8)   implicit none
     9)   
    10)   private 
    11) 
    12) #include "petsc/finclude/petscsys.h"
    13)   
    14) #include "petsc/finclude/petscvec.h"
    15) #include "petsc/finclude/petscvec.h90"
    16) #include "petsc/finclude/petscmat.h"
    17) #include "petsc/finclude/petscmat.h90"
    18) #include "petsc/finclude/petscsnes.h"
    19) #include "petsc/finclude/petscviewer.h"
    20) #include "petsc/finclude/petsclog.h"
    21) 
    22) 
    23) ! Cutoff parameters
    24)   PetscReal, parameter :: eps       = 1.D-8
    25)   PetscReal, parameter :: floweps   = 1.D-24
    26)   PetscReal, parameter :: perturbation_tolerance = 1.d-6
    27)   
    28)   public RichardsAccumulation, &
    29)          RichardsAccumDerivative, &
    30)          RichardsFlux, &
    31)          RichardsFluxDerivative, &
    32)          RichardsBCFlux, &
    33)          RichardsBCFluxDerivative
    34)   
    35) contains
    36) 
    37) ! ************************************************************************** !
    38) 
    39) subroutine RichardsAccumDerivative(rich_auxvar,global_auxvar, &
    40)                                    material_auxvar, &
    41)                                    option, &
    42)                                    characteristic_curves, &
    43)                                    J)
    44)   ! 
    45)   ! Computes derivatives of the accumulation
    46)   ! term for the Jacobian
    47)   ! 
    48)   ! Author: Glenn Hammond
    49)   ! Date: 12/13/07
    50)   ! 
    51) 
    52)   use Option_module
    53)   use Characteristic_Curves_module
    54)   use Material_Aux_class, only : material_auxvar_type, &
    55)                                  soil_compressibility_index, &
    56)                                  MaterialAuxVarInit, &
    57)                                  MaterialAuxVarCopy, &
    58)                                  MaterialAuxVarStrip, &
    59)                                  MaterialCompressSoil
    60)   
    61)   implicit none
    62) 
    63)   type(richards_auxvar_type) :: rich_auxvar
    64)   type(global_auxvar_type) :: global_auxvar
    65)   class(material_auxvar_type) :: material_auxvar
    66)   type(option_type) :: option
    67)   type(characteristic_curves_type) :: characteristic_curves
    68)   PetscReal :: J(option%nflowdof,option%nflowdof)
    69)      
    70)   PetscInt :: ispec 
    71)   PetscReal :: vol_over_dt
    72)   PetscReal :: por
    73)   PetscReal :: tempreal
    74)   PetscReal :: derivative
    75)   PetscReal :: compressed_porosity, dcompressed_porosity_dp
    76) 
    77)   PetscInt :: iphase, ideriv
    78)   type(richards_auxvar_type) :: rich_auxvar_pert
    79)   type(global_auxvar_type) :: global_auxvar_pert
    80)   ! leave as type
    81)   type(material_auxvar_type) :: material_auxvar_pert
    82)   PetscReal :: x(1), x_pert(1), pert, res(1), res_pert(1), J_pert(1,1)
    83) 
    84)   vol_over_dt = material_auxvar%volume/option%flow_dt
    85)   por = material_auxvar%porosity
    86)   
    87)   derivative = 0.d0
    88)   if (soil_compressibility_index > 0) then
    89)     tempreal = global_auxvar%sat(1)*global_auxvar%den(1)
    90)     call MaterialCompressSoil(material_auxvar,global_auxvar%pres(1), &
    91)                                  compressed_porosity,dcompressed_porosity_dp)
    92)     por = compressed_porosity
    93)     derivative = derivative + dcompressed_porosity_dp * tempreal
    94)   endif
    95) 
    96)   ! accumulation term units = dkmol/dp
    97)   derivative = derivative + (global_auxvar%sat(1)*rich_auxvar%dden_dp + &
    98)                              rich_auxvar%dsat_dp*global_auxvar%den(1)) * por
    99) 
   100)   J(1,1) = derivative * vol_over_dt
   101) 
   102)   if (option%flow%numerical_derivatives) then
   103)     call GlobalAuxVarInit(global_auxvar_pert,option)  
   104)     call MaterialAuxVarInit(material_auxvar_pert,option)  
   105)     call RichardsAuxVarCopy(rich_auxvar,rich_auxvar_pert,option)
   106)     call GlobalAuxVarCopy(global_auxvar,global_auxvar_pert,option)
   107)     call MaterialAuxVarCopy(material_auxvar,material_auxvar_pert,option)
   108)     x(1) = global_auxvar%pres(1)
   109)     call RichardsAccumulation(rich_auxvar,global_auxvar,material_auxvar, &
   110)                               option,res)
   111)     ideriv = 1
   112)     pert = max(dabs(x(ideriv)*perturbation_tolerance),0.1d0)
   113)     x_pert = x
   114)     if (x_pert(ideriv) < option%reference_pressure) pert = -1.d0*pert
   115)     x_pert(ideriv) = x_pert(ideriv) + pert
   116)     
   117)     call RichardsAuxVarCompute(x_pert(1),rich_auxvar_pert,global_auxvar_pert, &
   118)                                material_auxvar_pert, &
   119)                                characteristic_curves, &
   120)                                option)
   121)     call RichardsAccumulation(rich_auxvar_pert,global_auxvar_pert, &
   122)                               material_auxvar_pert, &
   123)                               option,res_pert)
   124)     J_pert(1,1) = (res_pert(1)-res(1))/pert
   125)     J = J_pert
   126)     call GlobalAuxVarStrip(global_auxvar_pert)  
   127)     call MaterialAuxVarStrip(material_auxvar_pert)  
   128)   endif
   129)    
   130) end subroutine RichardsAccumDerivative
   131) 
   132) ! ************************************************************************** !
   133) 
   134) subroutine RichardsAccumulation(rich_auxvar,global_auxvar, &
   135)                                 material_auxvar, &
   136)                                 option,Res)
   137)   ! 
   138)   ! Computes the non-fixed portion of the accumulation
   139)   ! term for the residual
   140)   ! 
   141)   ! Author: Glenn Hammond
   142)   ! Date: 12/13/07
   143)   ! 
   144) 
   145)   use Option_module
   146)   use Material_Aux_class, only : material_auxvar_type, &
   147)                                  soil_compressibility_index, &
   148)                                  MaterialCompressSoil
   149)   
   150)   implicit none
   151) 
   152)   type(richards_auxvar_type) :: rich_auxvar
   153)   type(global_auxvar_type) :: global_auxvar
   154)   class(material_auxvar_type) :: material_auxvar
   155)   type(option_type) :: option
   156)   PetscReal :: Res(1:option%nflowdof)
   157)   
   158)   PetscReal :: por, por1, vol_over_dt
   159)   PetscReal :: compressed_porosity, dcompressed_porosity_dp
   160)        
   161)   vol_over_dt = material_auxvar%volume/option%flow_dt
   162)   por = material_auxvar%porosity
   163)   
   164)   if (soil_compressibility_index > 0) then
   165)     call MaterialCompressSoil(material_auxvar,global_auxvar%pres(1), &
   166)                               compressed_porosity,dcompressed_porosity_dp)
   167)     material_auxvar%porosity = compressed_porosity
   168)     por = compressed_porosity
   169)   endif
   170)     
   171)   ! accumulation term units = kmol/s
   172)   Res(1) = global_auxvar%sat(1) * global_auxvar%den(1) * por * &
   173)            vol_over_dt
   174)   
   175) end subroutine RichardsAccumulation
   176) 
   177) ! ************************************************************************** !
   178) 
   179) subroutine RichardsFluxDerivative(rich_auxvar_up,global_auxvar_up, &
   180)                                   material_auxvar_up,sir_up, & 
   181)                                   rich_auxvar_dn,global_auxvar_dn, &
   182)                                   material_auxvar_dn,sir_dn, &
   183)                                   area, dist, &
   184)                                   option, &
   185)                                   characteristic_curves_up, &
   186)                                   characteristic_curves_dn, &
   187)                                   Jup,Jdn)
   188)   ! 
   189)   ! Computes the derivatives of the internal flux terms
   190)   ! for the Jacobian
   191)   ! 
   192)   ! Author: Glenn Hammond
   193)   ! Date: 12/13/07
   194)   ! 
   195)   use Option_module 
   196)   use Characteristic_Curves_module
   197)   use Material_Aux_class
   198)   use Connection_module
   199)   
   200)   implicit none
   201)   
   202)   type(richards_auxvar_type) :: rich_auxvar_up, rich_auxvar_dn
   203)   type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
   204)   class(material_auxvar_type) :: material_auxvar_up, material_auxvar_dn
   205)   type(option_type) :: option
   206)   PetscReal :: sir_up, sir_dn
   207)   PetscReal :: v_darcy, area, dist(-1:3)
   208)   type(characteristic_curves_type) :: characteristic_curves_up
   209)   type(characteristic_curves_type) :: characteristic_curves_dn
   210)   PetscReal :: Jup(option%nflowdof,option%nflowdof)
   211)   PetscReal :: Jdn(option%nflowdof,option%nflowdof)
   212)      
   213)   PetscReal :: q
   214)   PetscReal :: ukvr,Dq
   215)   PetscReal :: dd_up, dd_dn, perm_up, perm_dn
   216)   PetscReal :: dist_gravity
   217)   PetscReal :: upweight,density_ave,cond,gravity,dphi
   218)   
   219)   PetscReal :: dden_ave_dp_up, dden_ave_dp_dn
   220)   PetscReal :: dgravity_dden_up, dgravity_dden_dn
   221)   PetscReal :: dphi_dp_up, dphi_dp_dn
   222)   PetscReal :: dukvr_dp_up, dukvr_dp_dn
   223) !  PetscReal :: dukvr_x_dp_up, dukvr_x_dp_dn
   224) !  PetscReal :: dukvr_y_dp_up, dukvr_y_dp_dn
   225) !  PetscReal :: dukvr_z_dp_up, dukvr_z_dp_dn
   226)   PetscReal :: dq_dp_up, dq_dp_dn
   227) 
   228)   PetscInt :: iphase, ideriv
   229)   type(richards_auxvar_type) :: rich_auxvar_pert_up, rich_auxvar_pert_dn
   230)   type(global_auxvar_type) :: global_auxvar_pert_up, global_auxvar_pert_dn
   231)   ! leave as type
   232)   type(material_auxvar_type) :: material_auxvar_pert_up, material_auxvar_pert_dn
   233)   PetscReal :: x_up(1), x_dn(1), x_pert_up(1), x_pert_dn(1), pert_up, pert_dn, &
   234)             res(1), res_pert_up(1), res_pert_dn(1), J_pert_up(1,1), J_pert_dn(1,1)
   235)   
   236)   v_darcy = 0.D0  
   237)   ukvr = 0.d0
   238)   
   239)   Jup = 0.d0
   240)   Jdn = 0.d0 
   241)   
   242)   dden_ave_dp_up = 0.d0
   243)   dden_ave_dp_dn = 0.d0
   244)   dgravity_dden_up = 0.d0
   245)   dgravity_dden_dn = 0.d0
   246)   dphi_dp_up = 0.d0
   247)   dphi_dp_dn = 0.d0
   248)   dukvr_dp_up = 0.d0
   249)   dukvr_dp_dn = 0.d0
   250)   dq_dp_up = 0.d0
   251)   dq_dp_dn = 0.d0
   252)   
   253)   call ConnectionCalculateDistances(dist,option%gravity,dd_up,dd_dn, &
   254)                                     dist_gravity,upweight)
   255)   call material_auxvar_up%PermeabilityTensorToScalar(dist,perm_up)
   256)   call material_auxvar_dn%PermeabilityTensorToScalar(dist,perm_dn)
   257)   
   258)   Dq = (perm_up * perm_dn)/(dd_up*perm_dn + dd_dn*perm_up)
   259)   
   260) ! Flow term
   261)   if (global_auxvar_up%sat(1) > sir_up .or. global_auxvar_dn%sat(1) > sir_dn) then
   262)     if (global_auxvar_up%sat(1) <eps) then 
   263)       upweight=0.d0
   264)     else if (global_auxvar_dn%sat(1) <eps) then 
   265)       upweight=1.d0
   266)     endif
   267)     density_ave = upweight*global_auxvar_up%den(1)+ &
   268)                   (1.D0-upweight)*global_auxvar_dn%den(1)
   269)     dden_ave_dp_up = upweight*rich_auxvar_up%dden_dp
   270)     dden_ave_dp_dn = (1.D0-upweight)*rich_auxvar_dn%dden_dp
   271) 
   272)     gravity = (upweight*global_auxvar_up%den(1) + &
   273)               (1.D0-upweight)*global_auxvar_dn%den(1)) &
   274)               * FMWH2O * dist_gravity
   275)     dgravity_dden_up = upweight*FMWH2O*dist_gravity
   276)     dgravity_dden_dn = (1.d0-upweight)*FMWH2O*dist_gravity
   277) 
   278)     dphi = global_auxvar_up%pres(1) - global_auxvar_dn%pres(1)  + gravity
   279)     dphi_dp_up = 1.d0 + dgravity_dden_up*rich_auxvar_up%dden_dp
   280)     dphi_dp_dn = -1.d0 + dgravity_dden_dn*rich_auxvar_dn%dden_dp
   281) 
   282)     if (dphi>=0.D0) then
   283) #ifdef USE_ANISOTROPIC_MOBILITY
   284)       if (dabs(dist(1))==1) then
   285)         ukvr = rich_auxvar_up%kvr_x
   286)         dukvr_dp_up = rich_auxvar_up%dkvr_x_dp
   287)       else if (dabs(dist(2))==1) then
   288)         ukvr = rich_auxvar_up%kvr_y
   289)         dukvr_dp_up = rich_auxvar_up%dkvr_y_dp
   290)       else if (dabs(dist(3))==1) then
   291)         ukvr = rich_auxvar_up%kvr_z
   292)         dukvr_dp_up = rich_auxvar_up%dkvr_z_dp
   293)       end if
   294) #else
   295)       ukvr = rich_auxvar_up%kvr
   296)       dukvr_dp_up = rich_auxvar_up%dkvr_dp
   297) #endif
   298)     else
   299) #ifdef USE_ANISOTROPIC_MOBILITY    
   300)       if (dabs(dist(1))==1) then
   301)         ukvr = rich_auxvar_dn%kvr_x
   302)         dukvr_dp_dn = rich_auxvar_dn%dkvr_x_dp
   303)       else if (dabs(dist(2))==1) then
   304)         ukvr = rich_auxvar_dn%kvr_y
   305)         dukvr_dp_dn = rich_auxvar_dn%dkvr_y_dp
   306)       else if (dabs(dist(3))==1) then
   307)         ukvr = rich_auxvar_dn%kvr_z
   308)         dukvr_dp_dn = rich_auxvar_dn%dkvr_z_dp
   309)       end if
   310) #else
   311)       ukvr = rich_auxvar_dn%kvr
   312)       dukvr_dp_dn = rich_auxvar_dn%dkvr_dp
   313) #endif
   314)     endif      
   315)    
   316)     if (ukvr>floweps) then
   317)       v_darcy= Dq * ukvr * dphi
   318)    
   319)       q = v_darcy * area
   320)       dq_dp_up = Dq*(dukvr_dp_up*dphi+ukvr*dphi_dp_up)*area
   321)       dq_dp_dn = Dq*(dukvr_dp_dn*dphi+ukvr*dphi_dp_dn)*area
   322)       
   323)       Jup(1,1) = (dq_dp_up*density_ave+q*dden_ave_dp_up)
   324)       Jdn(1,1) = (dq_dp_dn*density_ave+q*dden_ave_dp_dn)
   325)     endif
   326)   endif
   327) 
   328)  ! note: Res is the flux contribution, for node up J = J + Jup
   329)  !                                              dn J = J - Jdn  
   330) 
   331)   if (option%flow%numerical_derivatives) then
   332)     call GlobalAuxVarInit(global_auxvar_pert_up,option)
   333)     call GlobalAuxVarInit(global_auxvar_pert_dn,option)  
   334)     call MaterialAuxVarInit(material_auxvar_pert_up,option)
   335)     call MaterialAuxVarInit(material_auxvar_pert_dn,option)  
   336)     call RichardsAuxVarCopy(rich_auxvar_up,rich_auxvar_pert_up,option)
   337)     call RichardsAuxVarCopy(rich_auxvar_dn,rich_auxvar_pert_dn,option)
   338)     call GlobalAuxVarCopy(global_auxvar_up,global_auxvar_pert_up,option)
   339)     call GlobalAuxVarCopy(global_auxvar_dn,global_auxvar_pert_dn,option)
   340)     call MaterialAuxVarCopy(material_auxvar_up,material_auxvar_pert_up,option)
   341)     call MaterialAuxVarCopy(material_auxvar_dn,material_auxvar_pert_dn,option)
   342)     x_up(1) = global_auxvar_up%pres(1)
   343)     x_dn(1) = global_auxvar_dn%pres(1)
   344)     call RichardsFlux(rich_auxvar_up,global_auxvar_up,material_auxvar_up, &
   345)                       sir_up, &
   346)                       rich_auxvar_dn,global_auxvar_dn,material_auxvar_dn, &
   347)                       sir_dn, &
   348)                       area, dist, &
   349)                       option,v_darcy,res)
   350)     ideriv = 1
   351) !    pert_up = x_up(ideriv)*perturbation_tolerance
   352)     pert_up = max(dabs(x_up(ideriv)*perturbation_tolerance),0.1d0)
   353)     if (x_up(ideriv) < option%reference_pressure) pert_up = -1.d0*pert_up
   354) !    pert_dn = x_dn(ideriv)*perturbation_tolerance
   355)     pert_dn = max(dabs(x_dn(ideriv)*perturbation_tolerance),0.1d0)
   356)     if (x_dn(ideriv) < option%reference_pressure) pert_dn = -1.d0*pert_dn
   357)     x_pert_up = x_up
   358)     x_pert_dn = x_dn
   359)     x_pert_up(ideriv) = x_pert_up(ideriv) + pert_up
   360)     x_pert_dn(ideriv) = x_pert_dn(ideriv) + pert_dn
   361)     call RichardsAuxVarCompute(x_pert_up(1),rich_auxvar_pert_up, &
   362)                                global_auxvar_pert_up, &
   363)                                material_auxvar_pert_up, &
   364)                                characteristic_curves_up, &
   365)                                option)
   366)     call RichardsAuxVarCompute(x_pert_dn(1),rich_auxvar_pert_dn, &
   367)                                global_auxvar_pert_dn, &
   368)                                material_auxvar_pert_dn, &
   369)                                characteristic_curves_dn, &
   370)                                option)
   371)     call RichardsFlux(rich_auxvar_pert_up,global_auxvar_pert_up, &
   372)                       material_auxvar_pert_up,sir_up, &
   373)                       rich_auxvar_dn,global_auxvar_dn, &
   374)                       material_auxvar_dn,sir_dn, &
   375)                       area, dist, &
   376)                       option,v_darcy,res_pert_up)
   377)     call RichardsFlux(rich_auxvar_up,global_auxvar_up, &
   378)                       material_auxvar_up,sir_up, &
   379)                       rich_auxvar_pert_dn,global_auxvar_pert_dn, &
   380)                       material_auxvar_pert_dn,sir_dn, &
   381)                       area, dist, &
   382)                       option,v_darcy,res_pert_dn)
   383)     J_pert_up(1,ideriv) = (res_pert_up(1)-res(1))/pert_up
   384)     J_pert_dn(1,ideriv) = (res_pert_dn(1)-res(1))/pert_dn
   385)     Jup = J_pert_up
   386)     Jdn = J_pert_dn
   387)     call GlobalAuxVarStrip(global_auxvar_pert_up)
   388)     call GlobalAuxVarStrip(global_auxvar_pert_dn)    
   389)     call MaterialAuxVarStrip(material_auxvar_pert_up)
   390)     call MaterialAuxVarStrip(material_auxvar_pert_dn)    
   391)   endif
   392) 
   393) end subroutine RichardsFluxDerivative
   394) 
   395) ! ************************************************************************** !
   396) 
   397) subroutine RichardsFlux(rich_auxvar_up,global_auxvar_up, &
   398)                         material_auxvar_up,sir_up, &
   399)                         rich_auxvar_dn,global_auxvar_dn, &
   400)                         material_auxvar_dn,sir_dn, &
   401)                         area, dist, &
   402)                         option,v_darcy,Res)
   403)   ! 
   404)   ! Computes the internal flux terms for the residual
   405)   ! 
   406)   ! Author: Glenn Hammond
   407)   ! Date: 12/13/07
   408)   ! 
   409)   use Option_module
   410)   use Material_Aux_class
   411)   use Connection_module
   412)   
   413)   implicit none
   414)   
   415)   type(richards_auxvar_type) :: rich_auxvar_up, rich_auxvar_dn
   416)   type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
   417)   class(material_auxvar_type) :: material_auxvar_up, material_auxvar_dn
   418)   type(option_type) :: option
   419)   PetscReal :: sir_up, sir_dn
   420)   PetscReal :: v_darcy, area, dist(-1:3)
   421)   PetscReal :: Res(1:option%nflowdof) 
   422)      
   423)   PetscInt :: ispec
   424)   PetscReal :: dist_gravity  ! distance along gravity vector
   425)   PetscReal :: dd_up, dd_dn, perm_up, perm_dn
   426)   PetscReal :: fluxm, q
   427)   PetscReal :: ukvr,Dq
   428)   PetscReal :: upweight, density_ave, cond, gravity, dphi
   429)   
   430)   fluxm = 0.d0
   431)   v_darcy = 0.D0  
   432)   ukvr = 0.d0
   433)   
   434)   call ConnectionCalculateDistances(dist,option%gravity,dd_up,dd_dn, &
   435)                                     dist_gravity,upweight)
   436)   call material_auxvar_up%PermeabilityTensorToScalar(dist,perm_up)
   437)   call material_auxvar_dn%PermeabilityTensorToScalar(dist,perm_dn)
   438) 
   439)   Dq = (perm_up * perm_dn)/(dd_up*perm_dn + dd_dn*perm_up)
   440)   
   441) ! Flow term
   442)   if (global_auxvar_up%sat(1) > sir_up .or. global_auxvar_dn%sat(1) > sir_dn) then
   443)     if (global_auxvar_up%sat(1) <eps) then 
   444)       upweight=0.d0
   445)     else if (global_auxvar_dn%sat(1) <eps) then 
   446)       upweight=1.d0
   447)     endif
   448)     
   449) 
   450)     density_ave = upweight*global_auxvar_up%den(1)+ &
   451)                   (1.D0-upweight)*global_auxvar_dn%den(1)
   452) 
   453)     gravity = (upweight*global_auxvar_up%den(1) + &
   454)               (1.D0-upweight)*global_auxvar_dn%den(1)) &
   455)               * FMWH2O * dist_gravity
   456) 
   457)     dphi = global_auxvar_up%pres(1) - global_auxvar_dn%pres(1)  + gravity
   458) 
   459)     if (dphi>=0.D0) then
   460) #ifdef USE_ANISOTROPIC_MOBILITY       
   461)       if (dabs(dabs(dist(1))-1) < 1e-6) then
   462)         ukvr = rich_auxvar_up%kvr_x
   463)       else if (dabs(dabs(norm(2))-1) < 1e-6) then
   464)         ukvr = rich_auxvar_up%kvr_y
   465)       else if (dabs(dabs(norm(3))-1) < 1e-6) then
   466)         ukvr = rich_auxvar_up%kvr_z
   467)       end if
   468) #else
   469)       ukvr = rich_auxvar_up%kvr
   470) #endif
   471)     else
   472) #ifdef USE_ANISOTROPIC_MOBILITY       
   473)       if (dabs(dabs(norm(1))-1) < 1e-6) then
   474)         ukvr = rich_auxvar_dn%kvr_x
   475)       else if (dabs(dabs(norm(2))-1) < 1e-6) then
   476)         ukvr = rich_auxvar_dn%kvr_y
   477)       else if (dabs(dabs(norm(3))-1) < 1e-6) then
   478)         ukvr = rich_auxvar_dn%kvr_z
   479)       end if
   480) #else
   481)       ukvr = rich_auxvar_dn%kvr
   482) #endif
   483)     endif      
   484) 
   485)     if (ukvr>floweps) then
   486)       v_darcy= Dq * ukvr * dphi
   487)    
   488)       q = v_darcy * area
   489) 
   490)       fluxm = q*density_ave       
   491)     endif
   492)   endif 
   493) 
   494)   Res(1) = fluxm
   495)  ! note: Res is the flux contribution, for node 1 R = R + Res_FL
   496)  !                                              2 R = R - Res_FL  
   497) 
   498) end subroutine RichardsFlux
   499) 
   500) ! ************************************************************************** !
   501) 
   502) subroutine RichardsBCFluxDerivative(ibndtype,auxvars, &
   503)                                     rich_auxvar_up,global_auxvar_up, &
   504)                                     rich_auxvar_dn,global_auxvar_dn, &
   505)                                     material_auxvar_dn, &
   506)                                     sir_dn, &
   507)                                     area,dist,option, &
   508)                                     characteristic_curves_dn, &
   509)                                     Jdn)
   510)   ! 
   511)   ! Computes the derivatives of the boundary flux
   512)   ! terms for the Jacobian
   513)   ! 
   514)   ! Author: Glenn Hammond
   515)   ! Date: 12/13/07
   516)   ! 
   517)   use Option_module
   518)   use Characteristic_Curves_module
   519)   use Material_Aux_class
   520)   use EOS_Water_module
   521)   use Utility_module
   522)  
   523)   implicit none
   524)   
   525)   PetscInt :: ibndtype(:)
   526)   type(richards_auxvar_type) :: rich_auxvar_up, rich_auxvar_dn
   527)   type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
   528)   class(material_auxvar_type) :: material_auxvar_dn
   529)   type(option_type) :: option
   530)   PetscReal :: sir_dn
   531)   PetscReal :: auxvars(:) ! from aux_real_var array in boundary condition
   532)   PetscReal :: area
   533)   ! dist(-1) = fraction_upwind
   534)   ! dist(0) = magnitude
   535)   ! dist(1:3) = unit vector
   536)   ! dist(0)*dist(1:3) = vector
   537)   PetscReal :: dist(-1:3)
   538)   type(characteristic_curves_type) :: characteristic_curves_dn
   539)   PetscReal :: Jdn(option%nflowdof,option%nflowdof)
   540)   
   541)   PetscReal :: dist_gravity  ! distance along gravity vector
   542)   PetscReal :: perm_dn
   543)   
   544)   PetscReal :: v_darcy
   545)   PetscReal :: q,density_ave
   546)   PetscReal :: ukvr,diffdp,Dq
   547)   PetscReal :: upweight,cond,gravity,dphi
   548) 
   549)   PetscReal :: dden_ave_dp_dn
   550)   PetscReal :: dgravity_dden_dn
   551)   PetscReal :: dphi_dp_dn
   552)   PetscReal :: dukvr_dp_dn
   553)   PetscReal :: dq_dp_dn
   554)   PetscInt :: pressure_bc_type
   555)   PetscReal :: dphi_x_dp_dn,dphi_y_dp_dn,dphi_z_dp_dn
   556)   PetscReal :: dHdn_x_dp_dn, dHdn_y_dp_dn, dHdn_z_dp_dn
   557) 
   558)   PetscInt :: iphase, ideriv
   559)   type(richards_auxvar_type) :: rich_auxvar_pert_dn, rich_auxvar_pert_up
   560)   type(global_auxvar_type) :: global_auxvar_pert_dn, global_auxvar_pert_up
   561)   class(material_auxvar_type), allocatable :: material_auxvar_pert_dn, &
   562)                                               material_auxvar_pert_up
   563)   PetscReal :: perturbation
   564)   PetscReal :: x_dn(1), x_up(1), x_pert_dn(1), x_pert_up(1), pert_dn, res(1), &
   565)             res_pert_dn(1), J_pert_dn(1,1)
   566)   PetscReal :: rho
   567)   PetscReal :: dum1
   568)   PetscReal :: dq_lin,dP_lin
   569)   PetscReal :: q_approx, dq_approx
   570)   PetscErrorCode :: ierr
   571) 
   572)   v_darcy = 0.d0
   573)   ukvr = 0.d0
   574)   density_ave = 0.d0
   575)   q = 0.d0
   576) 
   577)   Jdn = 0.d0 
   578)   
   579)   dden_ave_dp_dn = 0.d0
   580)   dgravity_dden_dn = 0.d0
   581)   dphi_dp_dn = 0.d0
   582)   dukvr_dp_dn = 0.d0
   583)   dq_dp_dn = 0.d0
   584) 
   585)   call material_auxvar_dn%PermeabilityTensorToScalar(dist,perm_dn)
   586)   
   587)   ! Flow
   588)   pressure_bc_type = ibndtype(RICHARDS_PRESSURE_DOF)
   589)   select case(pressure_bc_type)
   590)     ! figure out the direction of flow
   591)     case(DIRICHLET_BC,HYDROSTATIC_BC,SEEPAGE_BC,CONDUCTANCE_BC,HET_SURF_SEEPAGE_BC, &
   592)          HET_DIRICHLET)
   593) 
   594)       ! dist(0) = scalar - magnitude of distance
   595)       ! gravity = vector(3)
   596)       ! dist(1:3) = vector(3) - unit vector
   597)       dist_gravity = dist(0) * dot_product(option%gravity,dist(1:3))
   598) 
   599)       if (pressure_bc_type == CONDUCTANCE_BC) then
   600)         Dq = auxvars(RICHARDS_CONDUCTANCE_DOF)
   601)       else
   602)         Dq = perm_dn / dist(0)
   603)       endif
   604)       ! Flow term
   605)       if (global_auxvar_up%sat(1) > sir_dn .or. global_auxvar_dn%sat(1) > sir_dn) then
   606)         upweight=1.D0
   607)         if (global_auxvar_up%sat(1) < eps) then 
   608)           upweight=0.d0
   609)         else if (global_auxvar_dn%sat(1) < eps) then 
   610)           upweight=1.d0
   611)         endif
   612)         
   613)         density_ave = upweight*global_auxvar_up%den(1) + (1.D0-upweight)*global_auxvar_dn%den(1)
   614)         dden_ave_dp_dn = (1.D0-upweight)*rich_auxvar_dn%dden_dp
   615) 
   616)         gravity = (upweight*global_auxvar_up%den(1) + &
   617)                   (1.D0-upweight)*global_auxvar_dn%den(1)) &
   618)                   * FMWH2O * dist_gravity
   619)         dgravity_dden_dn = (1.d0-upweight)*FMWH2O*dist_gravity
   620) 
   621)         dphi = global_auxvar_up%pres(1) - global_auxvar_dn%pres(1) + gravity
   622)         dphi_dp_dn = -1.d0 + dgravity_dden_dn*rich_auxvar_dn%dden_dp
   623) 
   624)         if (pressure_bc_type == SEEPAGE_BC .or. &
   625)             pressure_bc_type == CONDUCTANCE_BC .or. &
   626)             pressure_bc_type == HET_SURF_SEEPAGE_BC) then
   627)               ! flow in         ! boundary cell is <= pref
   628)           if (dphi > 0.d0 .and. global_auxvar_up%pres(1)-option%reference_pressure < eps) then
   629)             dphi = 0.d0
   630)             dphi_dp_dn = 0.d0
   631)           endif
   632)         endif
   633)         
   634)         if (dphi>=0.D0) then
   635) #ifdef USE_ANISOTROPIC_MOBILITY  
   636)           if (dabs(dabs(dist(1))-1) < 1e-6) then
   637)             ukvr = rich_auxvar_up%kvr_x
   638)           else if (dabs(dabs(dist(2))-1) < 1e-6) then
   639)             ukvr = rich_auxvar_up%kvr_y
   640)           else if (dabs(dabs(dist(3))-1) < 1e-6) then
   641)             ukvr = rich_auxvar_up%kvr_z
   642)           end if
   643) #else
   644)           ukvr = rich_auxvar_up%kvr
   645) #endif
   646)         else
   647) #ifdef USE_ANISOTROPIC_MOBILITY
   648)           if (dabs(dabs(dist(1))-1) < 1e-6) then
   649)             ukvr = rich_auxvar_dn%kvr_x
   650)             dukvr_dp_dn = rich_auxvar_dn%dkvr_x_dp
   651)           else if (dabs(dabs(dist(2))-1) < 1e-6) then
   652)             ukvr = rich_auxvar_dn%kvr_y
   653)             dukvr_dp_dn = rich_auxvar_dn%dkvr_y_dp
   654)           else if (dabs(dabs(dist(3))-1) < 1e-6) then
   655)             ukvr = rich_auxvar_dn%kvr_z
   656)             dukvr_dp_dn = rich_auxvar_dn%dkvr_z_dp
   657)           end if
   658) #else
   659)           ukvr = rich_auxvar_dn%kvr
   660)           dukvr_dp_dn = rich_auxvar_dn%dkvr_dp
   661) #endif
   662)         endif      
   663)      
   664)         if (ukvr*Dq>floweps) then
   665)           v_darcy = Dq * ukvr * dphi
   666)           q = v_darcy * area
   667)           dq_dp_dn = Dq*(dukvr_dp_dn*dphi + ukvr*dphi_dp_dn)*area
   668) 
   669)           ! If running with surface-flow model, ensure (darcy_velocity*dt) does
   670)           ! not exceed depth of standing water.
   671)           if (option%surf_flow_on) then
   672)           if (rich_auxvar_dn%vars_for_sflow(11) == 0.d0) then
   673)             if (pressure_bc_type == HET_SURF_SEEPAGE_BC .and. option%surf_flow_on) then
   674)               call EOSWaterdensity(option%reference_temperature, &
   675)                                    option%reference_pressure,rho,dum1,ierr)
   676) 
   677)               if (global_auxvar_dn%pres(1) <= rich_auxvar_dn%vars_for_sflow(1)) then
   678) 
   679)                 ! Linear approximation
   680)                 call Interpolate(rich_auxvar_dn%vars_for_sflow(8), &
   681)                                  rich_auxvar_dn%vars_for_sflow(7), &
   682)                                  global_auxvar_dn%pres(1), &
   683)                                  rich_auxvar_dn%vars_for_sflow(10), &
   684)                                  rich_auxvar_dn%vars_for_sflow(9), &
   685)                                  q_approx)
   686)                 v_darcy = q_approx/area
   687)                 q       = q_approx
   688) 
   689)                 dP_lin = rich_auxvar_dn%vars_for_sflow(8) - &
   690)                          rich_auxvar_dn%vars_for_sflow(7)
   691)                 dq_lin = rich_auxvar_dn%vars_for_sflow(10) - &
   692)                          rich_auxvar_dn%vars_for_sflow(9)
   693)                 dq_dp_dn = dq_lin/dP_lin
   694) 
   695)               else
   696)                 if (global_auxvar_dn%pres(1) <= rich_auxvar_dn%vars_for_sflow(2)) then
   697) 
   698)                   ! Cubic approximation
   699)                   call CubicPolynomialEvaluate(rich_auxvar_dn%vars_for_sflow(3:6), &
   700)                                                global_auxvar_dn%pres(1) - option%reference_pressure, &
   701)                                                q_approx, dq_approx)
   702)                   v_darcy = q_approx/area
   703)                   q = q_approx
   704)                   dq_dp_dn = dq_approx
   705)                 endif
   706)               endif
   707)             endif
   708)           endif
   709)           endif
   710) 
   711)         endif
   712) 
   713)       endif
   714) 
   715)     case(NEUMANN_BC)
   716)       if (dabs(auxvars(RICHARDS_PRESSURE_DOF)) > floweps) then
   717)         v_darcy = auxvars(RICHARDS_PRESSURE_DOF)
   718)         if (v_darcy > 0.d0) then 
   719)           density_ave = global_auxvar_up%den(1)
   720)         else 
   721)           density_ave = global_auxvar_dn%den(1)
   722)           dden_ave_dp_dn = rich_auxvar_dn%dden_dp
   723)         endif 
   724)         q = v_darcy * area
   725)       endif
   726) 
   727)     case(UNIT_GRADIENT_BC)
   728) 
   729)       if (global_auxvar_dn%sat(1) > sir_dn) then
   730) 
   731)         dphi = dot_product(option%gravity,dist(1:3))*global_auxvar_dn%den(1)*FMWH2O      
   732)         density_ave = global_auxvar_dn%den(1)
   733) 
   734)         dphi_dp_dn = dot_product(option%gravity,dist(1:3))*rich_auxvar_dn%dden_dp*FMWH2O
   735)         dden_ave_dp_dn = rich_auxvar_dn%dden_dp
   736) 
   737)         ! since boundary auxvar is meaningless (no pressure specified there), only use cell
   738) #ifdef USE_ANISOTROPIC_MOBILITY
   739)         if (dabs(dabs(dist(1))-1) < 1e-6) then
   740)           ukvr = rich_auxvar_dn%kvr_x
   741)           dukvr_dp_dn = rich_auxvar_dn%dkvr_x_dp
   742)         else if (dabs(dabs(dist(2))-1) < 1e-6) then
   743)           ukvr = rich_auxvar_dn%kvr_y
   744)           dukvr_dp_dn = rich_auxvar_dn%dkvr_y_dp
   745)         else if (dabs(dabs(dist(3))-1) < 1e-6) then
   746)           ukvr = rich_auxvar_dn%kvr_z
   747)           dukvr_dp_dn = rich_auxvar_dn%dkvr_z_dp
   748)         end if
   749) #else
   750)         ukvr = rich_auxvar_dn%kvr
   751)         dukvr_dp_dn = rich_auxvar_dn%dkvr_dp
   752) #endif
   753)      
   754)         if (ukvr*perm_dn>floweps) then
   755)           v_darcy = perm_dn * ukvr * dphi
   756)           q = v_darcy*area
   757)           dq_dp_dn = perm_dn*(dukvr_dp_dn*dphi+ukvr*dphi_dp_dn)*area
   758)         endif 
   759)      
   760)       endif
   761) 
   762)   end select
   763) 
   764)   Jdn(1,1) = (dq_dp_dn*density_ave+q*dden_ave_dp_dn)
   765) 
   766)   if (option%flow%numerical_derivatives) then
   767)     call GlobalAuxVarInit(global_auxvar_pert_up,option)
   768)     call GlobalAuxVarInit(global_auxvar_pert_dn,option)  
   769)     allocate(material_auxvar_pert_up,material_auxvar_pert_dn)
   770)     call MaterialAuxVarInit(material_auxvar_pert_up,option)  
   771)     call MaterialAuxVarInit(material_auxvar_pert_dn,option)  
   772)     call RichardsAuxVarCopy(rich_auxvar_up,rich_auxvar_pert_up,option)
   773)     call RichardsAuxVarCopy(rich_auxvar_dn,rich_auxvar_pert_dn,option)
   774)     call GlobalAuxVarCopy(global_auxvar_up,global_auxvar_pert_up,option)
   775)     call GlobalAuxVarCopy(global_auxvar_dn,global_auxvar_pert_dn,option)
   776)     call MaterialAuxVarCopy(material_auxvar_dn,material_auxvar_pert_up, &
   777)                             option)
   778)     call MaterialAuxVarCopy(material_auxvar_dn,material_auxvar_pert_dn, &
   779)                             option)
   780)     x_up(1) = global_auxvar_up%pres(1)
   781)     x_dn(1) = global_auxvar_dn%pres(1)
   782)     ideriv = 1
   783)     if (ibndtype(ideriv) == ZERO_GRADIENT_BC) then
   784)       x_up(ideriv) = x_dn(ideriv)
   785)     endif
   786)     call RichardsBCFlux(ibndtype,auxvars, &
   787)                         rich_auxvar_up,global_auxvar_up, &
   788)                         rich_auxvar_dn,global_auxvar_dn, &
   789)                         material_auxvar_dn, &
   790)                         sir_dn, &
   791)                         area,dist,option,v_darcy,res)
   792)     if (pressure_bc_type == ZERO_GRADIENT_BC) then
   793)       x_pert_up = x_up
   794)     endif
   795)     ideriv = 1
   796) !    pert_dn = x_dn(ideriv)*perturbation_tolerance    
   797)     pert_dn = max(dabs(x_dn(ideriv)*perturbation_tolerance),0.1d0)
   798)     if (x_dn(ideriv) < option%reference_pressure) pert_dn = -1.d0*pert_dn
   799)     x_pert_dn = x_dn
   800)     x_pert_dn(ideriv) = x_pert_dn(ideriv) + pert_dn
   801)     x_pert_up = x_up
   802)     if (ibndtype(ideriv) == ZERO_GRADIENT_BC) then
   803)       x_pert_up(ideriv) = x_pert_dn(ideriv)
   804)     endif   
   805)     call RichardsAuxVarCompute(x_pert_dn(1),rich_auxvar_pert_dn, &
   806)                                global_auxvar_pert_dn, &
   807)                                material_auxvar_pert_dn, &
   808)                                characteristic_curves_dn, &
   809)                                option)
   810)     call RichardsAuxVarCompute(x_pert_up(1),rich_auxvar_pert_up, &
   811)                                global_auxvar_pert_up, &
   812)                                material_auxvar_pert_up, &
   813)                                characteristic_curves_dn, &
   814)                                option)
   815)     call RichardsBCFlux(ibndtype,auxvars, &
   816)                         rich_auxvar_pert_up,global_auxvar_pert_up, &
   817)                         rich_auxvar_pert_dn,global_auxvar_pert_dn, &
   818)                         material_auxvar_pert_dn, &
   819)                         sir_dn, &
   820)                         area,dist,option,v_darcy,res_pert_dn)
   821)     J_pert_dn(1,ideriv) = (res_pert_dn(1)-res(1))/pert_dn
   822)     Jdn = J_pert_dn
   823)     call GlobalAuxVarStrip(global_auxvar_pert_up)
   824)     call GlobalAuxVarStrip(global_auxvar_pert_dn)   
   825)     call MaterialAuxVarStrip(material_auxvar_pert_up)
   826)     call MaterialAuxVarStrip(material_auxvar_pert_dn)
   827)     deallocate(material_auxvar_pert_up,material_auxvar_pert_dn)
   828)   endif
   829) 
   830) end subroutine RichardsBCFluxDerivative
   831) 
   832) ! ************************************************************************** !
   833) 
   834) subroutine RichardsBCFlux(ibndtype,auxvars, &
   835)                           rich_auxvar_up, global_auxvar_up, &
   836)                           rich_auxvar_dn, global_auxvar_dn, &
   837)                           material_auxvar_dn, &
   838)                           sir_dn, &
   839)                           area, dist, option,v_darcy,Res)
   840)   ! 
   841)   ! Computes the  boundary flux terms for the residual
   842)   ! 
   843)   ! Author: Glenn Hammond
   844)   ! Date: 12/13/07
   845)   ! 
   846)   use Option_module
   847)   use Material_Aux_class
   848)   use EOS_Water_module
   849)   use Utility_module
   850)  
   851)   implicit none
   852)   
   853)   PetscInt :: ibndtype(:)
   854)   type(richards_auxvar_type) :: rich_auxvar_up, rich_auxvar_dn
   855)   type(global_auxvar_type) :: global_auxvar_up, global_auxvar_dn
   856)   class(material_auxvar_type) :: material_auxvar_dn
   857)   type(option_type) :: option
   858)   PetscReal :: sir_dn
   859)   PetscReal :: auxvars(:) ! from aux_real_var array
   860)   PetscReal :: v_darcy, area
   861)   ! dist(-1) = fraction_upwind - not applicable here
   862)   ! dist(0) = magnitude
   863)   ! dist(1:3) = unit vector
   864)   ! dist(0)*dist(1:3) = vector
   865)   PetscReal :: dist(-1:3)
   866)   PetscReal :: Res(1:option%nflowdof) 
   867)   
   868)   PetscReal :: dist_gravity  ! distance along gravity vector
   869)           
   870)   PetscInt :: ispec
   871)   PetscReal :: perm_dn
   872)   PetscReal :: fluxm,q,density_ave
   873)   PetscReal :: ukvr,diffdp,Dq
   874)   PetscReal :: upweight,cond,gravity,dphi
   875)   PetscInt :: pressure_bc_type
   876)   PetscReal :: dphi_x,dphi_y,dphi_z
   877)   PetscReal :: rho
   878)   PetscReal :: dum1
   879)   PetscReal :: q_approx, dq_approx
   880)   PetscErrorCode :: ierr
   881)   
   882)   fluxm = 0.d0
   883)   v_darcy = 0.d0
   884)   density_ave = 0.d0
   885)   q = 0.d0
   886)   ukvr = 0.d0
   887) 
   888)   call material_auxvar_dn%PermeabilityTensorToScalar(dist,perm_dn)  
   889) 
   890)   ! Flow  
   891)   pressure_bc_type = ibndtype(RICHARDS_PRESSURE_DOF)
   892)   select case(pressure_bc_type)
   893)     ! figure out the direction of flow
   894)     case(DIRICHLET_BC,HYDROSTATIC_BC,SEEPAGE_BC,CONDUCTANCE_BC,HET_SURF_SEEPAGE_BC, &
   895)          HET_DIRICHLET)
   896) 
   897)       ! dist(0) = scalar - magnitude of distance
   898)       ! gravity = vector(3)
   899)       ! dist(1:3) = vector(3) - unit vector
   900)       dist_gravity = dist(0) * dot_product(option%gravity,dist(1:3))
   901) 
   902)       if (pressure_bc_type == CONDUCTANCE_BC) then
   903)         Dq = auxvars(RICHARDS_CONDUCTANCE_DOF)
   904)       else
   905)         Dq = perm_dn / dist(0)
   906)       endif
   907)       ! Flow term
   908)       if (global_auxvar_up%sat(1) > sir_dn .or. global_auxvar_dn%sat(1) > sir_dn) then
   909)         upweight=1.D0
   910)         if (global_auxvar_up%sat(1) < eps) then 
   911)           upweight=0.d0
   912)         else if (global_auxvar_dn%sat(1) < eps) then 
   913)           upweight=1.d0
   914)         endif
   915)         density_ave = upweight*global_auxvar_up%den(1)+(1.D0-upweight)*global_auxvar_dn%den(1)
   916)    
   917)         gravity = (upweight*global_auxvar_up%den(1) + &
   918)                   (1.D0-upweight)*global_auxvar_dn%den(1)) &
   919)                   * FMWH2O * dist_gravity
   920)        
   921)         dphi = global_auxvar_up%pres(1) - global_auxvar_dn%pres(1) + gravity
   922) 
   923)         if (pressure_bc_type == SEEPAGE_BC .or. &
   924)             pressure_bc_type == CONDUCTANCE_BC .or. &
   925)             pressure_bc_type == HET_SURF_SEEPAGE_BC) then
   926)               ! flow in         ! boundary cell is <= pref
   927)           if (dphi > 0.d0 .and. global_auxvar_up%pres(1)-option%reference_pressure < eps) then
   928)             dphi = 0.d0
   929)           endif
   930)         endif
   931)    
   932)        if (dphi>=0.D0) then
   933) #ifdef USE_ANISOTROPIC_MOBILITY       
   934)          if (dabs(dabs(dist(1))-1) < 1e-6) then
   935)            ukvr = rich_auxvar_up%kvr_x
   936)          else if (dabs(dabs(dist(2))-1) < 1e-6) then
   937)            ukvr = rich_auxvar_up%kvr_y
   938)          else if (dabs(dabs(dist(3))-1) < 1e-6) then
   939)            ukvr = rich_auxvar_up%kvr_z
   940)          end if
   941) #else
   942)          ukvr = rich_auxvar_up%kvr
   943) #endif
   944)        else
   945) #ifdef USE_ANISOTROPIC_MOBILITY
   946)          if (dabs(dabs(dist(1))-1) < 1e-6) then
   947)            ukvr = rich_auxvar_dn%kvr_x
   948)          else if (dabs(dabs(dist(2))-1) < 1e-6) then
   949)            ukvr = rich_auxvar_dn%kvr_y
   950)          else if (dabs(dabs(dist(3))-1) < 1e-6) then
   951)            ukvr = rich_auxvar_dn%kvr_z
   952)          end if
   953) #else
   954)          ukvr = rich_auxvar_dn%kvr
   955) #endif
   956)        endif
   957)         
   958)        if (ukvr*Dq>floweps) then
   959)         v_darcy = Dq * ukvr * dphi
   960) 
   961)         ! If running with surface-flow model, ensure (darcy_velocity*dt) does
   962)         ! not exceed depth of standing water.
   963)         if (pressure_bc_type == HET_SURF_SEEPAGE_BC .and. option%surf_flow_on) then
   964)           call EOSWaterdensity(option%reference_temperature, &
   965)                                option%reference_pressure,rho,dum1,ierr)
   966) 
   967)           if (rich_auxvar_dn%vars_for_sflow(11) == 0.d0) then
   968)             if (global_auxvar_dn%pres(1) <= rich_auxvar_dn%vars_for_sflow(1)) then
   969) 
   970)               if (rich_auxvar_dn%vars_for_sflow(7) == -99999.d0) then
   971)                 call printErrMsg(option,'Coeffs for linear approx for darcy flux not set')
   972)               endif
   973) 
   974)               ! Linear approximation
   975)               call Interpolate(rich_auxvar_dn%vars_for_sflow(8), &
   976)                                rich_auxvar_dn%vars_for_sflow(7), &
   977)                                global_auxvar_dn%pres(1), &
   978)                                rich_auxvar_dn%vars_for_sflow(2), &
   979)                                rich_auxvar_dn%vars_for_sflow(1), &
   980)                                q_approx)
   981)               v_darcy = q_approx/area
   982) 
   983)             else if (global_auxvar_dn%pres(1) <= rich_auxvar_dn%vars_for_sflow(2)) then
   984) 
   985)               if (rich_auxvar_dn%vars_for_sflow(3) == -99999.d0) then
   986)                 call printErrMsg(option,'Coeffs for cubic approx for darcy flux not set')
   987)               endif
   988) 
   989)               ! Cubic approximation
   990)               call CubicPolynomialEvaluate(rich_auxvar_dn%vars_for_sflow(3:6), &
   991)                                            global_auxvar_dn%pres(1) - option%reference_pressure, &
   992)                                            q_approx, dq_approx)
   993)               v_darcy = q_approx/area
   994)             endif
   995)           endif
   996) 
   997)         endif
   998)        endif
   999)       endif 
  1000) 
  1001)     case(NEUMANN_BC)
  1002)       if (dabs(auxvars(RICHARDS_PRESSURE_DOF)) > floweps) then
  1003)         v_darcy = auxvars(RICHARDS_PRESSURE_DOF)
  1004) 
  1005)         if (v_darcy > 0.d0) then 
  1006)           density_ave = global_auxvar_up%den(1)
  1007)         else 
  1008)           density_ave = global_auxvar_dn%den(1)
  1009)         endif 
  1010)       endif
  1011) 
  1012)     case(UNIT_GRADIENT_BC)
  1013) 
  1014)       dphi = dot_product(option%gravity,dist(1:3))*global_auxvar_dn%den(1)*FMWH2O
  1015)       density_ave = global_auxvar_dn%den(1)
  1016) 
  1017)       ! since boundary auxvar is meaningless (no pressure specified there), only use cell
  1018) #ifdef USE_ANISOTROPIC_MOBILITY
  1019)       if (dabs(dabs(dist(1))-1) < 1e-6) then
  1020)         ukvr = rich_auxvar_dn%kvr_x
  1021)       else if (dabs(dabs(dist(2))-1) < 1e-6) then
  1022)         ukvr = rich_auxvar_dn%kvr_y
  1023)       else if (dabs(dabs(dist(3))-1) < 1e-6) then
  1024)         ukvr = rich_auxvar_dn%kvr_z
  1025)       end if
  1026) #else
  1027)       ukvr = rich_auxvar_dn%kvr
  1028) #endif
  1029)      
  1030)       if (ukvr*perm_dn>floweps) then
  1031)         v_darcy = perm_dn * ukvr * dphi
  1032)       endif 
  1033) 
  1034)   end select
  1035) 
  1036)   q = v_darcy * area
  1037) 
  1038)   fluxm = q*density_ave
  1039) 
  1040)   Res(1)=fluxm
  1041) 
  1042) end subroutine RichardsBCFlux
  1043) 
  1044) end module Richards_Common_module

generated by
Intel(R) C++/Fortran Compiler code-coverage tool
Web-Page Owner: Nobody