MODFLOW 6  version 6.5.0.dev2
MODFLOW 6 Code Documentation
gwegweexchangemodule Module Reference

This module contains the GweGweExchangeModule Module. More...

Data Types

type  gweexchangetype
 Derived type for GwtExchangeType. More...
 

Functions/Subroutines

subroutine, public gweexchange_create (filename, name, id, m1_id, m2_id, input_mempath)
 @ brief Create GWT GWT exchange More...
 
subroutine gwe_gwe_df (this)
 @ brief Define GWE GWE exchange More...
 
subroutine validate_exchange (this)
 validate exchange data after reading More...
 
subroutine gwe_gwe_ar (this)
 @ brief Allocate and read More...
 
subroutine gwe_gwe_rp (this)
 @ brief Read and prepare More...
 
subroutine gwe_gwe_ad (this)
 @ brief Advance More...
 
subroutine gwe_gwe_fc (this, kiter, matrix_sln, rhs_sln, inwtflag)
 @ brief Fill coefficients More...
 
subroutine gwe_gwe_bd (this, icnvg, isuppress_output, isolnid)
 @ brief Budget More...
 
subroutine gwe_gwe_bdsav (this)
 @ brief Budget save More...
 
subroutine gwe_gwe_bdsav_model (this, model)
 @ brief Budget save More...
 
subroutine gwe_gwe_ot (this)
 @ brief Output More...
 
subroutine source_options (this, iout)
 @ brief Source options More...
 
subroutine read_mvt (this, iout)
 @ brief Read mover More...
 
subroutine allocate_scalars (this)
 @ brief Allocate scalars More...
 
subroutine gwe_gwe_da (this)
 @ brief Deallocate More...
 
subroutine allocate_arrays (this)
 @ brief Allocate arrays More...
 
subroutine gwe_gwe_df_obs (this)
 @ brief Define observations More...
 
subroutine gwe_gwe_rp_obs (this)
 @ brief Read and prepare observations More...
 
subroutine gwe_gwe_fp (this)
 @ brief Final processing More...
 
logical(lgp) function gwe_gwe_connects_model (this, model)
 Return true when this exchange provides matrix coefficients for solving. More...
 
logical(lgp) function use_interface_model (this)
 Should interface model be used for this exchange. More...
 
subroutine gwe_gwe_save_simvals (this)
 @ brief Save simulated flow observations More...
 
subroutine gwe_gwe_process_obsid (obsrv, dis, inunitobs, iout)
 @ brief Obs ID processor More...
 
class(gweexchangetype) function, pointer, public castasgweexchange (obj)
 @ brief Cast polymorphic object as exchange More...
 
class(gweexchangetype) function, pointer, public getgweexchangefromlist (list, idx)
 @ brief Get exchange from list More...
 

Detailed Description

This module contains the code for connecting two GWE Models. The methods are based on the simple two point flux approximation with the option to use ghost nodes to improve accuracy. This exchange is used by GweGweConnection with the more sophisticated interface model coupling approach when XT3D is needed.

Function/Subroutine Documentation

◆ allocate_arrays()

subroutine gwegweexchangemodule::allocate_arrays ( class(gweexchangetype this)

Allocate arrays

Parameters
thisGweExchangeType

Definition at line 892 of file exg-gwegwe.f90.

893  ! -- modules
895  ! -- dummy
896  class(GweExchangeType) :: this !< GweExchangeType
897  ! -- local
898  character(len=LINELENGTH) :: text
899  integer(I4B) :: ntabcol, i
900  !
901  call this%DisConnExchangeType%allocate_arrays()
902  !
903  call mem_allocate(this%cond, this%nexg, 'COND', this%memoryPath)
904  call mem_allocate(this%simvals, this%nexg, 'SIMVALS', this%memoryPath)
905  !
906  ! -- Initialize
907  do i = 1, this%nexg
908  this%cond(i) = dnodata
909  end do
910  !
911  ! -- allocate and initialize the output table
912  if (this%iprflow /= 0) then
913  !
914  ! -- dimension table
915  ntabcol = 3
916  if (this%inamedbound > 0) then
917  ntabcol = ntabcol + 1
918  end if
919  !
920  ! -- initialize the output table objects
921  ! outouttab1
922  if (this%v_model1%is_local) then
923  call table_cr(this%outputtab1, this%name, ' ')
924  call this%outputtab1%table_df(this%nexg, ntabcol, this%gwemodel1%iout, &
925  transient=.true.)
926  text = 'NUMBER'
927  call this%outputtab1%initialize_column(text, 10, alignment=tabcenter)
928  text = 'CELLID'
929  call this%outputtab1%initialize_column(text, 20, alignment=tableft)
930  text = 'RATE'
931  call this%outputtab1%initialize_column(text, 15, alignment=tabcenter)
932  if (this%inamedbound > 0) then
933  text = 'NAME'
934  call this%outputtab1%initialize_column(text, 20, alignment=tableft)
935  end if
936  end if
937  ! outouttab2
938  if (this%v_model2%is_local) then
939  call table_cr(this%outputtab2, this%name, ' ')
940  call this%outputtab2%table_df(this%nexg, ntabcol, this%gwemodel2%iout, &
941  transient=.true.)
942  text = 'NUMBER'
943  call this%outputtab2%initialize_column(text, 10, alignment=tabcenter)
944  text = 'CELLID'
945  call this%outputtab2%initialize_column(text, 20, alignment=tableft)
946  text = 'RATE'
947  call this%outputtab2%initialize_column(text, 15, alignment=tabcenter)
948  if (this%inamedbound > 0) then
949  text = 'NAME'
950  call this%outputtab2%initialize_column(text, 20, alignment=tableft)
951  end if
952  end if
953  end if
954  !
955  ! -- Return
956  return
Here is the call graph for this function:

◆ allocate_scalars()

subroutine gwegweexchangemodule::allocate_scalars ( class(gweexchangetype this)

Allocate scalar variables

Parameters
thisGwtExchangeType

Definition at line 816 of file exg-gwegwe.f90.

817  ! -- modules
819  use constantsmodule, only: dzero
820  ! -- dummy
821  class(GweExchangeType) :: this !< GwtExchangeType
822  !
823  call this%DisConnExchangeType%allocate_scalars()
824  !
825  call mem_allocate(this%inewton, 'INEWTON', this%memoryPath)
826  call mem_allocate(this%inobs, 'INOBS', this%memoryPath)
827  call mem_allocate(this%iAdvScheme, 'IADVSCHEME', this%memoryPath)
828  this%inewton = 0
829  this%inobs = 0
830  this%iAdvScheme = 0
831  !
832  call mem_allocate(this%inmvt, 'INMVT', this%memoryPath)
833  this%inmvt = 0
834  !
835  ! -- Return
836  return
This module contains simulation constants.
Definition: Constants.f90:9
real(dp), parameter dzero
real constant zero
Definition: Constants.f90:64

◆ castasgweexchange()

class(gweexchangetype) function, pointer, public gwegweexchangemodule::castasgweexchange ( class(*), intent(inout), pointer  obj)

Cast polymorphic object as exchange

Definition at line 1205 of file exg-gwegwe.f90.

1206  implicit none
1207  ! -- dummy
1208  class(*), pointer, intent(inout) :: obj
1209  ! -- return
1210  class(GweExchangeType), pointer :: res
1211  !
1212  res => null()
1213  if (.not. associated(obj)) return
1214  !
1215  select type (obj)
1216  class is (gweexchangetype)
1217  res => obj
1218  end select
1219  !
1220  ! -- Return
1221  return
Here is the caller graph for this function:

◆ getgweexchangefromlist()

class(gweexchangetype) function, pointer, public gwegweexchangemodule::getgweexchangefromlist ( type(listtype), intent(inout)  list,
integer(i4b), intent(in)  idx 
)

Return an exchange from the list for specified index

Definition at line 1228 of file exg-gwegwe.f90.

1229  implicit none
1230  ! -- dummy
1231  type(ListType), intent(inout) :: list
1232  integer(I4B), intent(in) :: idx
1233  ! -- return
1234  class(GweExchangeType), pointer :: res
1235  ! -- local
1236  class(*), pointer :: obj
1237  !
1238  obj => list%GetItem(idx)
1239  res => castasgweexchange(obj)
1240  !
1241  ! -- Return
1242  return
Here is the call graph for this function:

◆ gwe_gwe_ad()

subroutine gwegweexchangemodule::gwe_gwe_ad ( class(gweexchangetype this)

Advance mover and obs

Parameters
thisGweExchangeType

Definition at line 358 of file exg-gwegwe.f90.

359  ! -- dummy
360  class(GweExchangeType) :: this !< GweExchangeType
361  !
362  ! -- Advance mover
363  !if(this%inmvt > 0) call this%mvt%mvt_ad()
364  !
365  ! -- Push simulated values to preceding time step
366  call this%obs%obs_ad()
367  !
368  ! -- Return
369  return

◆ gwe_gwe_ar()

subroutine gwegweexchangemodule::gwe_gwe_ar ( class(gweexchangetype this)
private

Allocated and read and calculate saturated conductance

Parameters
thisGwtExchangeType

Definition at line 317 of file exg-gwegwe.f90.

318  ! -- dummy
319  class(GweExchangeType) :: this !< GwtExchangeType
320  !
321  ! -- If mover is active, then call ar routine
322  if (this%inmvt > 0) call this%mvt%mvt_ar()
323  !
324  ! -- Observation AR
325  call this%obs%obs_ar()
326  !
327  ! -- Return
328  return

◆ gwe_gwe_bd()

subroutine gwegweexchangemodule::gwe_gwe_bd ( class(gweexchangetype this,
integer(i4b), intent(inout)  icnvg,
integer(i4b), intent(in)  isuppress_output,
integer(i4b), intent(in)  isolnid 
)
private

Accumulate budget terms

Parameters
thisGweExchangeType

Definition at line 395 of file exg-gwegwe.f90.

396  ! -- modules
398  use budgetmodule, only: rate_accumulator
399  ! -- dummy
400  class(GweExchangeType) :: this !< GweExchangeType
401  integer(I4B), intent(inout) :: icnvg
402  integer(I4B), intent(in) :: isuppress_output
403  integer(I4B), intent(in) :: isolnid
404  ! -- local
405  character(len=LENBUDTXT), dimension(1) :: budtxt
406  real(DP), dimension(2, 1) :: budterm
407  real(DP) :: ratin, ratout
408  !
409  ! -- initialize
410  budtxt(1) = ' FLOW-JA-FACE'
411  !
412  ! -- Calculate ratin/ratout and pass to model budgets
413  call rate_accumulator(this%simvals, ratin, ratout)
414  !
415  ! -- Add the budget terms to model 1
416  if (associated(this%gwemodel1)) then
417  budterm(1, 1) = ratin
418  budterm(2, 1) = ratout
419  call this%gwemodel1%model_bdentry(budterm, budtxt, this%name)
420  end if
421  !
422  ! -- Add the budget terms to model 2
423  if (associated(this%gwemodel2)) then
424  budterm(1, 1) = ratout
425  budterm(2, 1) = ratin
426  call this%gwemodel2%model_bdentry(budterm, budtxt, this%name)
427  end if
428  !
429  ! -- Call mvt bd routine
430  if (this%inmvt > 0) call this%mvt%mvt_bd(this%gwemodel1%x, this%gwemodel2%x)
431  !
432  ! -- Return
433  return
This module contains the BudgetModule.
Definition: Budget.f90:20
subroutine, public rate_accumulator(flow, rin, rout)
@ brief Rate accumulator subroutine
Definition: Budget.f90:664
integer(i4b), parameter lenpackagename
maximum length of the package name
Definition: Constants.f90:22
integer(i4b), parameter lenbudtxt
maximum length of a budget component names
Definition: Constants.f90:36
Here is the call graph for this function:

◆ gwe_gwe_bdsav()

subroutine gwegweexchangemodule::gwe_gwe_bdsav ( class(gweexchangetype this)

Output individual flows to listing file and binary budget files

Parameters
thisGweExchangeType

Definition at line 440 of file exg-gwegwe.f90.

441  ! -- dummy
442  class(GweExchangeType) :: this !< GweExchangeType
443  ! -- local
444  integer(I4B) :: icbcfl, ibudfl
445  !
446  ! -- budget for model1
447  if (associated(this%gwemodel1)) then
448  call this%gwe_gwe_bdsav_model(this%gwemodel1)
449  end if
450  !
451  ! -- budget for model2
452  if (associated(this%gwemodel2)) then
453  call this%gwe_gwe_bdsav_model(this%gwemodel2)
454  end if
455  !
456  ! -- Set icbcfl, ibudfl to zero so that flows will be printed and
457  ! saved, if the options were set in the MVT package
458  icbcfl = 1
459  ibudfl = 1
460  !
461  ! -- Call mvt bd routine
462  !cdl todo: if(this%inmvt > 0) call this%mvt%mvt_bdsav(icbcfl, ibudfl, isuppress_output)
463  !
464  ! -- Calculate and write simulated values for observations
465  if (this%inobs /= 0) then
466  call this%gwe_gwe_save_simvals()
467  end if
468  !
469  ! -- Return
470  return

◆ gwe_gwe_bdsav_model()

subroutine gwegweexchangemodule::gwe_gwe_bdsav_model ( class(gweexchangetype this,
class(gwemodeltype), pointer  model 
)
private

Output individual flows to listing file and binary budget files

Parameters
thisGwtExchangeType

Definition at line 477 of file exg-gwegwe.f90.

478  ! -- modules
480  use tdismodule, only: kstp, kper
481  ! -- dummy
482  class(GweExchangeType) :: this !< GwtExchangeType
483  class(GweModelType), pointer :: model
484  ! -- local
485  character(len=LENBOUNDNAME) :: bname
486  character(len=LENPACKAGENAME + 4) :: packname
487  character(len=LENBUDTXT), dimension(1) :: budtxt
488  type(TableType), pointer :: output_tab
489  class(VirtualModelType), pointer :: nbr_model
490  character(len=20) :: nodestr
491  integer(I4B) :: ntabrows
492  integer(I4B) :: nodeu
493  integer(I4B) :: i, n1, n2, n1u, n2u
494  integer(I4B) :: ibinun
495  real(DP) :: ratin, ratout, rrate
496  logical(LGP) :: is_for_model1
497  integer(I4B) :: isuppress_output
498  !
499  ! -- initialize local variables
500  isuppress_output = 0
501  budtxt(1) = ' FLOW-JA-FACE'
502  packname = 'EXG '//this%name
503  packname = adjustr(packname)
504  if (associated(model, this%gwemodel1)) then
505  output_tab => this%outputtab1
506  nbr_model => this%v_model2
507  is_for_model1 = .true.
508  else
509  output_tab => this%outputtab2
510  nbr_model => this%v_model1
511  is_for_model1 = .false.
512  end if
513  !
514  ! -- update output tables
515  if (this%iprflow /= 0) then
516  !
517  ! -- update titles
518  if (model%oc%oc_save('BUDGET')) then
519  call output_tab%set_title(packname)
520  end if
521  !
522  ! -- set table kstp and kper
523  call output_tab%set_kstpkper(kstp, kper)
524  !
525  ! -- update maxbound of tables
526  ntabrows = 0
527  do i = 1, this%nexg
528  n1 = this%nodem1(i)
529  n2 = this%nodem2(i)
530  !
531  ! -- If both cells are active then calculate flow rate
532  if (this%v_model1%ibound%get(n1) /= 0 .and. &
533  this%v_model2%ibound%get(n2) /= 0) then
534  ntabrows = ntabrows + 1
535  end if
536  end do
537  if (ntabrows > 0) then
538  call output_tab%set_maxbound(ntabrows)
539  end if
540  end if
541  !
542  ! -- Print and write budget terms for model 1
543  !
544  ! -- Set binary unit numbers for saving flows
545  if (this%ipakcb /= 0) then
546  ibinun = model%oc%oc_save_unit('BUDGET')
547  else
548  ibinun = 0
549  end if
550  !
551  ! -- If save budget flag is zero for this stress period, then
552  ! shut off saving
553  if (.not. model%oc%oc_save('BUDGET')) ibinun = 0
554  if (isuppress_output /= 0) then
555  ibinun = 0
556  end if
557  !
558  ! -- If cell-by-cell flows will be saved as a list, write header.
559  if (ibinun /= 0) then
560  call model%dis%record_srcdst_list_header(budtxt(1), &
561  model%name, &
562  this%name, &
563  nbr_model%name, &
564  this%name, &
565  this%naux, this%auxname, &
566  ibinun, this%nexg, &
567  model%iout)
568  end if
569  !
570  ! Initialize accumulators
571  ratin = dzero
572  ratout = dzero
573  !
574  ! -- Loop through all exchanges
575  do i = 1, this%nexg
576  !
577  ! -- Assign boundary name
578  if (this%inamedbound > 0) then
579  bname = this%boundname(i)
580  else
581  bname = ''
582  end if
583  !
584  ! -- Calculate the flow rate between n1 and n2
585  rrate = dzero
586  n1 = this%nodem1(i)
587  n2 = this%nodem2(i)
588  !
589  ! -- If both cells are active then calculate flow rate
590  if (this%v_model1%ibound%get(n1) /= 0 .and. &
591  this%v_model2%ibound%get(n2) /= 0) then
592  rrate = this%simvals(i)
593  !
594  ! -- Print the individual rates to model list files if requested
595  if (this%iprflow /= 0) then
596  if (model%oc%oc_save('BUDGET')) then
597  !
598  ! -- set nodestr and write outputtab table
599  if (is_for_model1) then
600  nodeu = model%dis%get_nodeuser(n1)
601  call model%dis%nodeu_to_string(nodeu, nodestr)
602  call output_tab%print_list_entry(i, trim(adjustl(nodestr)), &
603  rrate, bname)
604  else
605  nodeu = model%dis%get_nodeuser(n2)
606  call model%dis%nodeu_to_string(nodeu, nodestr)
607  call output_tab%print_list_entry(i, trim(adjustl(nodestr)), &
608  -rrate, bname)
609  end if
610  end if
611  end if
612  if (rrate < dzero) then
613  ratout = ratout - rrate
614  else
615  ratin = ratin + rrate
616  end if
617  end if
618  !
619  ! -- If saving cell-by-cell flows in list, write flow
620  n1u = this%v_model1%dis_get_nodeuser(n1)
621  n2u = this%v_model2%dis_get_nodeuser(n2)
622  if (ibinun /= 0) then
623  if (is_for_model1) then
624  call model%dis%record_mf6_list_entry( &
625  ibinun, n1u, n2u, rrate, this%naux, this%auxvar(:, i), &
626  .false., .false.)
627  else
628  call model%dis%record_mf6_list_entry( &
629  ibinun, n2u, n1u, -rrate, this%naux, this%auxvar(:, i), &
630  .false., .false.)
631  end if
632  end if
633  !
634  end do
635  !
636  ! -- Return
637  return
integer(i4b), pointer, public kstp
current time step number
Definition: tdis.f90:24
integer(i4b), pointer, public kper
current stress period number
Definition: tdis.f90:23

◆ gwe_gwe_connects_model()

logical(lgp) function gwegweexchangemodule::gwe_gwe_connects_model ( class(gweexchangetype this,
class(basemodeltype), intent(in), pointer  model 
)
private
Parameters
model
thisGweExchangeType
[in]modelthe model to which the exchange might hold a connection
Returns
true, when connected

Definition at line 1067 of file exg-gwegwe.f90.

1068  ! -- dummy
1069  class(GweExchangeType) :: this !< GweExchangeType
1070  class(BaseModelType), pointer, intent(in) :: model !< the model to which the exchange might hold a connection
1071  ! -- return
1072  logical(LGP) :: is_connected !< true, when connected
1073  !
1074  is_connected = .false.
1075  !
1076  ! only connected when model is GwtModelType of course
1077  select type (model)
1078  class is (gwemodeltype)
1079  if (associated(this%gwemodel1, model)) then
1080  is_connected = .true.
1081  else if (associated(this%gwemodel2, model)) then
1082  is_connected = .true.
1083  end if
1084  end select
1085  !
1086  ! -- Return
1087  return

◆ gwe_gwe_da()

subroutine gwegweexchangemodule::gwe_gwe_da ( class(gweexchangetype this)

Deallocate memory associated with this object

Parameters
thisGwtExchangeType

Definition at line 843 of file exg-gwegwe.f90.

844  ! -- modules
846  ! -- dummy
847  class(GweExchangeType) :: this !< GwtExchangeType
848  !
849  ! -- objects
850  if (this%inmvt > 0) then
851  call this%mvt%mvt_da()
852  deallocate (this%mvt)
853  end if
854  call this%obs%obs_da()
855  deallocate (this%obs)
856  !
857  ! -- arrays
858  call mem_deallocate(this%cond)
859  call mem_deallocate(this%simvals)
860  call mem_deallocate(this%gwfsimvals, 'GWFSIMVALS', this%memoryPath) ! linked memory
861  !
862  ! -- output table objects
863  if (associated(this%outputtab1)) then
864  call this%outputtab1%table_da()
865  deallocate (this%outputtab1)
866  nullify (this%outputtab1)
867  end if
868  if (associated(this%outputtab2)) then
869  call this%outputtab2%table_da()
870  deallocate (this%outputtab2)
871  nullify (this%outputtab2)
872  end if
873  !
874  ! -- scalars
875  deallocate (this%filename)
876  call mem_deallocate(this%inewton)
877  call mem_deallocate(this%inobs)
878  call mem_deallocate(this%iAdvScheme)
879  call mem_deallocate(this%inmvt)
880  !
881  ! -- deallocate base
882  call this%DisConnExchangeType%disconnex_da()
883  !
884  ! -- Return
885  return

◆ gwe_gwe_df()

subroutine gwegweexchangemodule::gwe_gwe_df ( class(gweexchangetype this)

Define GWE to GWE exchange object.

Parameters
thisGwtExchangeType

Definition at line 198 of file exg-gwegwe.f90.

199  ! -- modules
200  use simvariablesmodule, only: iout
202  use ghostnodemodule, only: gnc_cr
203  ! -- dummy
204  class(GweExchangeType) :: this !< GwtExchangeType
205  !
206  ! -- log the exchange
207  write (iout, '(/a,a)') ' Creating exchange: ', this%name
208  !
209  ! -- Ensure models are in same solution
210  if (associated(this%gwemodel1) .and. associated(this%gwemodel2)) then
211  if (this%gwemodel1%idsoln /= this%gwemodel2%idsoln) then
212  call store_error('Two models are connect in a GWE '// &
213  'exchange but they are in different solutions. '// &
214  'GWE models must be in same solution: '// &
215  trim(this%gwemodel1%name)//' '// &
216  trim(this%gwemodel2%name))
217  call store_error_filename(this%filename)
218  end if
219  end if
220  !
221  ! -- source options
222  call this%source_options(iout)
223  !
224  ! -- source dimensions
225  call this%source_dimensions(iout)
226  !
227  ! -- allocate arrays
228  call this%allocate_arrays()
229  !
230  ! -- source exchange data
231  call this%source_data(iout)
232  !
233  ! -- Read mover information
234  if (this%inmvt > 0) then
235  call this%read_mvt(iout)
236  call this%mvt%mvt_df(this%gwemodel1%dis)
237  end if
238  !
239  ! -- Store obs
240  call this%gwe_gwe_df_obs()
241  if (associated(this%gwemodel1)) then
242  call this%obs%obs_df(iout, this%name, 'GWE-GWE', this%gwemodel1%dis)
243  end if
244  !
245  ! -- validate
246  call this%validate_exchange()
247  !
248  ! -- Return
249  return
subroutine, public gnc_cr(gncobj, name_parent, inunit, iout)
Create new GNC exchange object.
Definition: GhostNode.f90:61
integer(i4b) function, public getunit()
Get a free unit number.
subroutine, public openfile(iu, iout, fname, ftype, fmtarg_opt, accarg_opt, filstat_opt, mode_opt)
Open a file.
Definition: InputOutput.f90:30
This module contains simulation variables.
Definition: SimVariables.f90:9
integer(i4b) iout
file unit number for simulation output
Here is the call graph for this function:

◆ gwe_gwe_df_obs()

subroutine gwegweexchangemodule::gwe_gwe_df_obs ( class(gweexchangetype this)

Define the observations associated with this object

Parameters
thisGweExchangeType

Definition at line 963 of file exg-gwegwe.f90.

964  ! -- dummy
965  class(GweExchangeType) :: this !< GweExchangeType
966  ! -- local
967  integer(I4B) :: indx
968  !
969  ! -- Store obs type and assign procedure pointer
970  ! for gwt-gwt observation type.
971  call this%obs%StoreObsType('flow-ja-face', .true., indx)
972  this%obs%obsData(indx)%ProcessIdPtr => gwe_gwe_process_obsid
973  !
974  ! -- Return
975  return
Here is the call graph for this function:

◆ gwe_gwe_fc()

subroutine gwegweexchangemodule::gwe_gwe_fc ( class(gweexchangetype this,
integer(i4b), intent(in)  kiter,
class(matrixbasetype), pointer  matrix_sln,
real(dp), dimension(:), intent(inout)  rhs_sln,
integer(i4b), intent(in), optional  inwtflag 
)
private

Calculate conductance and fill coefficient matrix

Parameters
thisGwtExchangeType

Definition at line 376 of file exg-gwegwe.f90.

377  ! -- dummy
378  class(GweExchangeType) :: this !< GwtExchangeType
379  integer(I4B), intent(in) :: kiter
380  class(MatrixBaseType), pointer :: matrix_sln
381  real(DP), dimension(:), intent(inout) :: rhs_sln
382  integer(I4B), optional, intent(in) :: inwtflag
383  !
384  ! -- Call mvt fc routine
385  if (this%inmvt > 0) call this%mvt%mvt_fc(this%gwemodel1%x, this%gwemodel2%x)
386  !
387  ! -- Return
388  return

◆ gwe_gwe_fp()

subroutine gwegweexchangemodule::gwe_gwe_fp ( class(gweexchangetype this)

Conduct any final processing

Parameters
thisGwtExchangeType

Definition at line 1056 of file exg-gwegwe.f90.

1057  ! -- dummy
1058  class(GweExchangeType) :: this !< GwtExchangeType
1059  !
1060  ! -- Return
1061  return

◆ gwe_gwe_ot()

subroutine gwegweexchangemodule::gwe_gwe_ot ( class(gweexchangetype this)

Write output

Parameters
thisGweExchangeType

Definition at line 644 of file exg-gwegwe.f90.

645  ! -- modules
646  use simvariablesmodule, only: iout
647  use constantsmodule, only: dzero
648  ! -- dummy
649  class(GweExchangeType) :: this !< GweExchangeType
650  ! -- local
651  integer(I4B) :: iexg, n1, n2
652  integer(I4B) :: ibudfl
653  real(DP) :: flow
654  character(len=LINELENGTH) :: node1str, node2str
655  ! -- format
656  character(len=*), parameter :: fmtheader = &
657  "(/1x, 'SUMMARY OF EXCHANGE RATES FOR EXCHANGE ', a, ' WITH ID ', i0, /, &
658  &2a16, 5a16, /, 112('-'))"
659  character(len=*), parameter :: fmtheader2 = &
660  "(/1x, 'SUMMARY OF EXCHANGE RATES FOR EXCHANGE ', a, ' WITH ID ', i0, /, &
661  &2a16, 4a16, /, 96('-'))"
662  character(len=*), parameter :: fmtdata = &
663  "(2a16, 5(1pg16.6))"
664  !
665  ! -- Call bdsave
666  call this%gwe_gwe_bdsav()
667  !
668  ! -- Write a table of exchanges
669  if (this%iprflow /= 0) then
670  write (iout, fmtheader2) trim(adjustl(this%name)), this%id, 'NODEM1', &
671  'NODEM2', 'COND', 'X_M1', 'X_M2', 'FLOW'
672  do iexg = 1, this%nexg
673  n1 = this%nodem1(iexg)
674  n2 = this%nodem2(iexg)
675  flow = this%simvals(iexg)
676  call this%v_model1%dis_noder_to_string(n1, node1str)
677  call this%v_model2%dis_noder_to_string(n2, node2str)
678  write (iout, fmtdata) trim(adjustl(node1str)), &
679  trim(adjustl(node2str)), &
680  this%cond(iexg), this%v_model1%x%get(n1), &
681  this%v_model2%x%get(n2), flow
682  end do
683  end if
684  !
685  !cdl Implement when MVT is ready
686  ! -- Mover budget output
687  ibudfl = 1
688  if (this%inmvt > 0) call this%mvt%mvt_ot_bdsummary(ibudfl)
689  !
690  ! -- OBS output
691  call this%obs%obs_ot()
692  !
693  ! -- Return
694  return

◆ gwe_gwe_process_obsid()

subroutine gwegweexchangemodule::gwe_gwe_process_obsid ( type(observetype), intent(inout)  obsrv,
class(disbasetype), intent(in)  dis,
integer(i4b), intent(in)  inunitobs,
integer(i4b), intent(in)  iout 
)

Process observations for this exchange

Definition at line 1163 of file exg-gwegwe.f90.

1164  ! -- modules
1165  use constantsmodule, only: linelength
1166  use inputoutputmodule, only: urword
1167  use observemodule, only: observetype
1168  use basedismodule, only: disbasetype
1169  ! -- dummy
1170  type(ObserveType), intent(inout) :: obsrv
1171  class(DisBaseType), intent(in) :: dis
1172  integer(I4B), intent(in) :: inunitobs
1173  integer(I4B), intent(in) :: iout
1174  ! -- local
1175  integer(I4B) :: n, iexg, istat
1176  integer(I4B) :: icol, istart, istop
1177  real(DP) :: r
1178  character(len=LINELENGTH) :: string
1179  !
1180  string = obsrv%IDstring
1181  icol = 1
1182  ! -- get exchange index
1183  call urword(string, icol, istart, istop, 1, n, r, iout, inunitobs)
1184  read (string(istart:istop), '(i10)', iostat=istat) iexg
1185  if (istat == 0) then
1186  obsrv%intPak1 = iexg
1187  else
1188  ! Integer can't be read from string; it's presumed to be an exchange
1189  ! boundary name (already converted to uppercase)
1190  obsrv%FeatureName = trim(adjustl(string))
1191  ! -- Observation may require summing rates from multiple exchange
1192  ! boundaries, so assign intPak1 as a value that indicates observation
1193  ! is for a named exchange boundary or group of exchange boundaries.
1194  obsrv%intPak1 = namedboundflag
1195  end if
1196  !
1197  ! -- Return
1198  return
integer(i4b), parameter linelength
maximum length of a standard line
Definition: Constants.f90:44
subroutine, public urword(line, icol, istart, istop, ncode, n, r, iout, in)
Extract a word from a string.
This module contains the derived types ObserveType and ObsDataType.
Definition: Observe.f90:15
Here is the call graph for this function:
Here is the caller graph for this function:

◆ gwe_gwe_rp()

subroutine gwegweexchangemodule::gwe_gwe_rp ( class(gweexchangetype this)
private

Read new data for mover and obs

Parameters
thisGweExchangeType

Definition at line 335 of file exg-gwegwe.f90.

336  ! -- modules
337  use tdismodule, only: readnewdata
338  ! -- dummy
339  class(GweExchangeType) :: this !< GweExchangeType
340  !
341  ! -- Check with TDIS on whether or not it is time to RP
342  if (.not. readnewdata) return
343  !
344  ! -- Read and prepare for mover
345  if (this%inmvt > 0) call this%mvt%mvt_rp()
346  !
347  ! -- Read and prepare for observations
348  call this%gwe_gwe_rp_obs()
349  !
350  ! -- Return
351  return
logical(lgp), pointer, public readnewdata
flag indicating time to read new data
Definition: tdis.f90:26

◆ gwe_gwe_rp_obs()

subroutine gwegweexchangemodule::gwe_gwe_rp_obs ( class(gweexchangetype this)
private

Handle observation exchanges exchange-boundary names.

Parameters
thisGwtExchangeType

Definition at line 982 of file exg-gwegwe.f90.

983  ! -- modules
984  use constantsmodule, only: dzero
985  ! -- dummy
986  class(GweExchangeType) :: this !< GwtExchangeType
987  ! -- local
988  integer(I4B) :: i
989  integer(I4B) :: j
990  class(ObserveType), pointer :: obsrv => null()
991  character(len=LENBOUNDNAME) :: bname
992  logical :: jfound
993  ! -- formats
994 10 format('Exchange "', a, '" for observation "', a, &
995  '" is invalid in package "', a, '"')
996 20 format('Exchange id "', i0, '" for observation "', a, &
997  '" is invalid in package "', a, '"')
998  !
999  do i = 1, this%obs%npakobs
1000  obsrv => this%obs%pakobs(i)%obsrv
1001  !
1002  ! -- indxbnds needs to be reset each stress period because
1003  ! list of boundaries can change each stress period.
1004  ! -- Not true for exchanges, but leave this in for now anyway.
1005  call obsrv%ResetObsIndex()
1006  obsrv%BndFound = .false.
1007  !
1008  bname = obsrv%FeatureName
1009  if (bname /= '') then
1010  ! -- Observation location(s) is(are) based on a boundary name.
1011  ! Iterate through all boundaries to identify and store
1012  ! corresponding index(indices) in bound array.
1013  jfound = .false.
1014  do j = 1, this%nexg
1015  if (this%boundname(j) == bname) then
1016  jfound = .true.
1017  obsrv%BndFound = .true.
1018  obsrv%CurrentTimeStepEndValue = dzero
1019  call obsrv%AddObsIndex(j)
1020  end if
1021  end do
1022  if (.not. jfound) then
1023  write (errmsg, 10) trim(bname), trim(obsrv%ObsTypeId), trim(this%name)
1024  call store_error(errmsg)
1025  end if
1026  else
1027  ! -- Observation location is a single exchange number
1028  if (obsrv%intPak1 <= this%nexg .and. obsrv%intPak1 > 0) then
1029  jfound = .true.
1030  obsrv%BndFound = .true.
1031  obsrv%CurrentTimeStepEndValue = dzero
1032  call obsrv%AddObsIndex(obsrv%intPak1)
1033  else
1034  jfound = .false.
1035  end if
1036  if (.not. jfound) then
1037  write (errmsg, 20) obsrv%intPak1, trim(obsrv%ObsTypeId), trim(this%name)
1038  call store_error(errmsg)
1039  end if
1040  end if
1041  end do
1042  !
1043  ! -- write summary of error messages
1044  if (count_errors() > 0) then
1045  call store_error_filename(this%obs%inputFilename)
1046  end if
1047  !
1048  ! -- Return
1049  return
Here is the call graph for this function:

◆ gwe_gwe_save_simvals()

subroutine gwegweexchangemodule::gwe_gwe_save_simvals ( class(gweexchangetype), intent(inout)  this)
private

Save the simulated flows for each exchange

Definition at line 1116 of file exg-gwegwe.f90.

1117  ! -- dummy
1118  use simvariablesmodule, only: errmsg
1119  use constantsmodule, only: dzero
1120  use observemodule, only: observetype
1121  class(GweExchangeType), intent(inout) :: this
1122  ! -- local
1123  integer(I4B) :: i
1124  integer(I4B) :: j
1125  integer(I4B) :: n1
1126  integer(I4B) :: n2
1127  integer(I4B) :: iexg
1128  real(DP) :: v
1129  type(ObserveType), pointer :: obsrv => null()
1130  !
1131  ! -- Write simulated values for all gwt-gwt observations
1132  if (this%obs%npakobs > 0) then
1133  call this%obs%obs_bd_clear()
1134  do i = 1, this%obs%npakobs
1135  obsrv => this%obs%pakobs(i)%obsrv
1136  do j = 1, obsrv%indxbnds_count
1137  iexg = obsrv%indxbnds(j)
1138  v = dzero
1139  select case (obsrv%ObsTypeId)
1140  case ('FLOW-JA-FACE')
1141  n1 = this%nodem1(iexg)
1142  n2 = this%nodem2(iexg)
1143  v = this%simvals(iexg)
1144  case default
1145  errmsg = 'Unrecognized observation type: '// &
1146  trim(obsrv%ObsTypeId)
1147  call store_error(errmsg)
1148  call store_error_filename(this%obs%inputFilename)
1149  end select
1150  call this%obs%SaveOneSimval(obsrv, v)
1151  end do
1152  end do
1153  end if
1154  !
1155  ! -- Return
1156  return
character(len=maxcharlen) errmsg
error message string
Here is the call graph for this function:

◆ gweexchange_create()

subroutine, public gwegweexchangemodule::gweexchange_create ( character(len=*), intent(in)  filename,
character(len=*)  name,
integer(i4b), intent(in)  id,
integer(i4b), intent(in)  m1_id,
integer(i4b), intent(in)  m2_id,
character(len=*), intent(in)  input_mempath 
)

Create a new GWT to GWT exchange object.

Parameters
[in]filenamefilename for reading
[in]idid for the exchange
namethe exchange name
[in]m1_idid for model 1
[in]m2_idid for model 2

Definition at line 110 of file exg-gwegwe.f90.

111  ! -- modules
112  use basemodelmodule, only: basemodeltype
113  use listsmodule, only: baseexchangelist
114  use obsmodule, only: obs_cr
116  ! -- dummy
117  character(len=*), intent(in) :: filename !< filename for reading
118  integer(I4B), intent(in) :: id !< id for the exchange
119  character(len=*) :: name !< the exchange name
120  integer(I4B), intent(in) :: m1_id !< id for model 1
121  integer(I4B), intent(in) :: m2_id !< id for model 2
122  character(len=*), intent(in) :: input_mempath
123  ! -- local
124  type(GweExchangeType), pointer :: exchange
125  class(BaseModelType), pointer :: mb
126  class(BaseExchangeType), pointer :: baseexchange
127  integer(I4B) :: m1_index, m2_index
128  !
129  ! -- Create a new exchange and add it to the baseexchangelist container
130  allocate (exchange)
131  baseexchange => exchange
132  call addbaseexchangetolist(baseexchangelist, baseexchange)
133  !
134  ! -- Assign id and name
135  exchange%id = id
136  exchange%name = name
137  exchange%memoryPath = create_mem_path(exchange%name)
138  exchange%input_mempath = input_mempath
139  !
140  ! -- allocate scalars and set defaults
141  call exchange%allocate_scalars()
142  exchange%filename = filename
143  exchange%typename = 'GWE-GWE'
144  exchange%iAdvScheme = 0
145  exchange%ixt3d = 1
146  !
147  ! -- set gwemodel1
148  m1_index = model_loc_idx(m1_id)
149  mb => getbasemodelfromlist(basemodellist, m1_index)
150  if (m1_index > 0) then
151  select type (mb)
152  type is (gwemodeltype)
153  exchange%model1 => mb
154  exchange%gwemodel1 => mb
155  end select
156  end if
157  exchange%v_model1 => get_virtual_model(m1_id)
158  !
159  ! -- set gwemodel2
160  m2_index = model_loc_idx(m2_id)
161  if (m2_index > 0) then
162  mb => getbasemodelfromlist(basemodellist, m2_index)
163  select type (mb)
164  type is (gwemodeltype)
165  exchange%model2 => mb
166  exchange%gwemodel2 => mb
167  end select
168  end if
169  exchange%v_model2 => get_virtual_model(m2_id)
170  !
171  ! -- Verify that gwt model1 is of the correct type
172  if (.not. associated(exchange%gwemodel1) .and. m1_index > 0) then
173  write (errmsg, '(3a)') 'Problem with GWE-GWE exchange ', &
174  trim(exchange%name), &
175  '. First specified GWE Model does not appear to be of the correct type.'
176  call store_error(errmsg, terminate=.true.)
177  end if
178  !
179  ! -- Verify that gwe model2 is of the correct type
180  if (.not. associated(exchange%gwemodel2) .and. m2_index > 0) then
181  write (errmsg, '(3a)') 'Problem with GWE-GWE exchange ', &
182  trim(exchange%name), &
183  '. Second specified GWE Model does not appear to be of the correct type.'
184  call store_error(errmsg, terminate=.true.)
185  end if
186  !
187  ! -- Create the obs package
188  call obs_cr(exchange%obs, exchange%inobs)
189  !
190  ! -- Return
191  return
type(listtype), public baseexchangelist
Definition: mf6lists.f90:25
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains the derived type ObsType.
Definition: Obs.f90:127
subroutine, public obs_cr(obs, inobs)
@ brief Create a new ObsType object
Definition: Obs.f90:225
Highest level model type. All models extend this parent type.
Definition: BaseModel.f90:13
Here is the call graph for this function:
Here is the caller graph for this function:

◆ read_mvt()

subroutine gwegweexchangemodule::read_mvt ( class(gweexchangetype this,
integer(i4b), intent(in)  iout 
)

Read and process movers

Parameters
thisGwtExchangeType

Definition at line 792 of file exg-gwegwe.f90.

793  ! -- modules
794  use tspmvtmodule, only: mvt_cr
795  ! -- dummy
796  class(GweExchangeType) :: this !< GwtExchangeType
797  integer(I4B), intent(in) :: iout
798  !
799  ! -- Create and initialize the mover object Here, fmi is set to the one
800  ! for gwtmodel1 so that a call to save flows has an associated dis
801  ! object.
802  call mvt_cr(this%mvt, this%name, this%inmvt, iout, this%gwemodel1%fmi, &
803  this%gwemodel1%eqnsclfac, this%gwemodel1%depvartype, &
804  gwfmodelname1=this%gwfmodelname1, &
805  gwfmodelname2=this%gwfmodelname2, &
806  fmi2=this%gwemodel2%fmi)
807  !
808  ! -- Return
809  return
subroutine, public mvt_cr(mvt, name_model, inunit, iout, fmi1, eqnsclfac, depvartype, gwfmodelname1, gwfmodelname2, fmi2)
Create a new mover transport object.
Definition: tsp-mvt.f90:75
Here is the call graph for this function:

◆ source_options()

subroutine gwegweexchangemodule::source_options ( class(gweexchangetype this,
integer(i4b), intent(in)  iout 
)

Source the options block

Parameters
thisGweExchangeType

Definition at line 701 of file exg-gwegwe.f90.

702  ! -- modules
703  use constantsmodule, only: lenvarname
709  ! -- dummy
710  class(GweExchangeType) :: this !< GweExchangeType
711  integer(I4B), intent(in) :: iout
712  ! -- local
713  type(ExgGwegweParamFoundType) :: found
714  character(len=LENVARNAME), dimension(3) :: adv_scheme = &
715  &[character(len=LENVARNAME) :: 'UPSTREAM', 'CENTRAL', 'TVD']
716  character(len=linelength) :: mvt_fname
717  !
718  ! -- update defaults with values sourced from input context
719  call mem_set_value(this%gwfmodelname1, 'GWFMODELNAME1', this%input_mempath, &
720  found%gwfmodelname1)
721  call mem_set_value(this%gwfmodelname2, 'GWFMODELNAME2', this%input_mempath, &
722  found%gwfmodelname2)
723  call mem_set_value(this%iAdvScheme, 'ADV_SCHEME', this%input_mempath, &
724  adv_scheme, found%adv_scheme)
725  call mem_set_value(this%ixt3d, 'CND_XT3D_OFF', this%input_mempath, &
726  found%cnd_xt3d_off)
727  call mem_set_value(this%ixt3d, 'CND_XT3D_RHS', this%input_mempath, &
728  found%cnd_xt3d_rhs)
729  !
730  write (iout, '(1x,a)') 'PROCESSING GWE-GWE EXCHANGE OPTIONS'
731  !
732  ! -- source base class options
733  call this%DisConnExchangeType%source_options(iout)
734  !
735  if (found%gwfmodelname1) then
736  write (iout, '(4x,a,a)') &
737  'GWFMODELNAME1 IS SET TO: ', trim(this%gwfmodelname1)
738  end if
739  !
740  if (found%gwfmodelname2) then
741  write (iout, '(4x,a,a)') &
742  'GWFMODELNAME2 IS SET TO: ', trim(this%gwfmodelname2)
743  end if
744  !
745  if (found%adv_scheme) then
746  ! -- count from 0
747  this%iAdvScheme = this%iAdvScheme - 1
748  write (iout, '(4x,a,a)') &
749  'ADVECTION SCHEME METHOD HAS BEEN SET TO: ', &
750  trim(adv_scheme(this%iAdvScheme + 1))
751  end if
752  !
753  if (found%cnd_xt3d_off .and. found%cnd_xt3d_rhs) then
754  errmsg = 'CND_XT3D_OFF and CND_XT3D_RHS cannot both be set as options.'
755  call store_error(errmsg)
756  call store_error_filename(this%filename)
757  else if (found%cnd_xt3d_off) then
758  this%ixt3d = 0
759  write (iout, '(4x,a)') 'XT3D FORMULATION HAS BEEN SHUT OFF.'
760  else if (found%cnd_xt3d_rhs) then
761  this%ixt3d = 2
762  write (iout, '(4x,a)') 'XT3D RIGHT-HAND SIDE FORMULATION IS SELECTED.'
763  end if
764  !
765  ! -- enforce 0 or 1 MVR6_FILENAME entries in option block
766  if (filein_fname(mvt_fname, 'MVE6_FILENAME', this%input_mempath, &
767  this%filename)) then
768  this%inmvt = getunit()
769  call openfile(this%inmvt, iout, mvt_fname, 'MVT')
770  write (iout, '(4x,a)') 'WATER MOVER ENERGY TRANSPORT &
771  &INFORMATION WILL BE READ FROM ', trim(mvt_fname)
772  end if
773  !
774  ! -- enforce 0 or 1 OBS6_FILENAME entries in option block
775  if (filein_fname(this%obs%inputFilename, 'OBS6_FILENAME', &
776  this%input_mempath, this%filename)) then
777  this%obs%active = .true.
778  this%obs%inUnitObs = getunit()
779  call openfile(this%obs%inUnitObs, iout, this%obs%inputFilename, 'OBS')
780  end if
781  !
782  write (iout, '(1x,a)') 'END OF GWE-GWE EXCHANGE OPTIONS'
783  !
784  ! -- return
785  return
integer(i4b), parameter lenvarname
maximum length of a variable name
Definition: Constants.f90:17
This module contains the SourceCommonModule.
Definition: SourceCommon.f90:7
logical(lgp) function, public filein_fname(filename, tagname, input_mempath, input_fname)
enforce and set a single input filename provided via FILEIN keyword
This class is used to store a single deferred-length character string. It was designed to work in an ...
Definition: CharString.f90:23
Here is the call graph for this function:

◆ use_interface_model()

logical(lgp) function gwegweexchangemodule::use_interface_model ( class(gweexchangetype this)
private

For now this always returns true, since we do not support a classic-style two-point flux approximation for GWT-GWT. If we ever add logic to support a simpler non-interface model flux calculation, then logic should be added here to set the return accordingly.

Parameters
thisGweExchangeType
Returns
true when interface model should be used

Definition at line 1098 of file exg-gwegwe.f90.

1099  ! -- dummy
1100  class(GweExchangeType) :: this !< GweExchangeType
1101  ! -- return
1102  logical(LGP) :: use_im !< true when interface model should be used
1103  !
1104  ! For now set use_im to .true. since the interface model approach
1105  ! must currently be used for any GWT-GWT exchange.
1106  use_im = .true.
1107  !
1108  ! -- Return
1109  return

◆ validate_exchange()

subroutine gwegweexchangemodule::validate_exchange ( class(gweexchangetype this)
Parameters
thisGweExchangeType

Definition at line 254 of file exg-gwegwe.f90.

255  ! -- dummy
256  class(GweExchangeType) :: this !< GweExchangeType
257  !
258 
259  ! Ensure gwfmodel names were entered
260  if (this%gwfmodelname1 == '') then
261  write (errmsg, '(3a)') 'GWE-GWE exchange ', trim(this%name), &
262  ' requires that GWFMODELNAME1 be entered in the &
263  &OPTIONS block.'
264  call store_error(errmsg)
265  end if
266  if (this%gwfmodelname2 == '') then
267  write (errmsg, '(3a)') 'GWE-GWE exchange ', trim(this%name), &
268  ' requires that GWFMODELNAME2 be entered in the &
269  &OPTIONS block.'
270  call store_error(errmsg)
271  end if
272  !
273  ! Periodic boundary condition in exchange don't allow XT3D (=interface model)
274  if (associated(this%model1, this%model2)) then
275  if (this%ixt3d > 0) then
276  write (errmsg, '(3a)') 'GWE-GWE exchange ', trim(this%name), &
277  ' is a periodic boundary condition which cannot'// &
278  ' be configured with XT3D'
279  call store_error(errmsg)
280  end if
281  end if
282  !
283  ! Check to see if dispersion is on in either model1 or model2.
284  ! If so, then ANGLDEGX must be provided as an auxiliary variable for this
285  ! GWE-GWE exchange (this%ianglex > 0).
286  if (associated(this%gwemodel1) .and. associated(this%gwemodel2)) then
287  if (this%gwemodel1%incnd /= 0 .or. this%gwemodel2%incnd /= 0) then
288  if (this%ianglex == 0) then
289  write (errmsg, '(3a)') 'GWE-GWE exchange ', trim(this%name), &
290  ' requires that ANGLDEGX be specified as an'// &
291  ' auxiliary variable because dispersion was '// &
292  'specified in one or both transport models.'
293  call store_error(errmsg)
294  end if
295  end if
296  end if
297  !
298  if (this%ixt3d > 0 .and. this%ianglex == 0) then
299  write (errmsg, '(3a)') 'GWE-GWE exchange ', trim(this%name), &
300  ' requires that ANGLDEGX be specified as an'// &
301  ' auxiliary variable because XT3D is enabled'
302  call store_error(errmsg)
303  end if
304  !
305  if (count_errors() > 0) then
306  call ustop()
307  end if
308  !
309  ! -- Return
310  return
Here is the call graph for this function: