15 character(len=LENMEMPATH) :: memorypath
16 integer(I4B),
pointer :: inunit => null()
17 integer(I4B),
pointer :: iout => null()
18 integer(I4B),
pointer :: inewton => null()
19 integer(I4B),
dimension(:),
pointer,
contiguous :: ibound => null()
20 integer(I4B),
dimension(:),
pointer,
contiguous :: iax => null()
21 integer(I4B),
dimension(:),
pointer,
contiguous :: jax => null()
22 integer(I4B),
dimension(:),
pointer,
contiguous :: idxglox => null()
23 integer(I4B),
dimension(:),
pointer,
contiguous :: ia_xt3d => null()
24 integer(I4B),
dimension(:),
pointer,
contiguous :: ja_xt3d => null()
25 integer(I4B),
pointer :: numextnbrs => null()
26 integer(I4B),
pointer :: ixt3d => null()
27 logical,
pointer :: nozee => null()
28 real(dp),
pointer :: vcthresh => null()
29 real(dp),
dimension(:, :),
pointer,
contiguous :: rmatck => null()
30 real(dp),
dimension(:),
pointer,
contiguous :: qsat => null()
31 integer(I4B),
pointer :: nbrmax => null()
32 real(dp),
dimension(:),
pointer,
contiguous :: amatpc => null()
33 real(dp),
dimension(:),
pointer,
contiguous :: amatpcx => null()
34 integer(I4B),
dimension(:),
pointer,
contiguous :: iallpc => null()
35 logical,
pointer :: lamatsaved => null()
38 real(dp),
dimension(:),
pointer,
contiguous :: k11 => null()
39 real(dp),
dimension(:),
pointer,
contiguous :: k22 => null()
40 real(dp),
dimension(:),
pointer,
contiguous :: k33 => null()
41 integer(I4B),
pointer :: ik22 => null()
42 integer(I4B),
pointer :: ik33 => null()
43 real(dp),
dimension(:),
pointer,
contiguous :: sat => null()
44 integer(I4B),
dimension(:),
pointer,
contiguous :: icelltype => null()
45 integer(I4B),
pointer :: iangle1 => null()
46 integer(I4B),
pointer :: iangle2 => null()
47 integer(I4B),
pointer :: iangle3 => null()
48 real(dp),
dimension(:),
pointer,
contiguous :: angle1 => null()
49 real(dp),
dimension(:),
pointer,
contiguous :: angle2 => null()
50 real(dp),
dimension(:),
pointer,
contiguous :: angle3 => null()
51 logical,
pointer :: ldispersion => null()
89 subroutine xt3d_cr(xt3dobj, name_model, inunit, iout, ldispopt)
92 character(len=*),
intent(in) :: name_model
93 integer(I4B),
intent(in) :: inunit
94 integer(I4B),
intent(in) :: iout
95 logical,
optional,
intent(in) :: ldispopt
104 call xt3dobj%allocate_scalars()
107 xt3dobj%inunit = inunit
109 if (
present(ldispopt)) xt3dobj%ldispersion = ldispopt
130 integer(I4B),
intent(in) :: moffset
134 integer(I4B) :: i, j, k, jj, kk, iglo, kglo, iadded
136 integer(I4B) :: ierror
139 if (this%ixt3d == 1)
then
143 call sparse_xt3d%init(this%dis%nodes, this%dis%nodes, nnz)
148 do i = 1, this%dis%nodes
151 do jj = this%dis%con%ia(i) + 1, this%dis%con%ia(i + 1) - 1
152 j = this%dis%con%ja(jj)
154 do kk = this%dis%con%ia(j) + 1, this%dis%con%ia(j + 1) - 1
155 k = this%dis%con%ja(kk)
157 call sparse_xt3d%addconnection(i, k, 1)
164 call mem_allocate(this%ia_xt3d, this%dis%nodes + 1,
'IA_XT3D', &
165 trim(this%memoryPath))
166 call mem_allocate(this%ja_xt3d, sparse_xt3d%nnz,
'JA_XT3D', &
167 trim(this%memoryPath))
168 call sparse_xt3d%filliaja(this%ia_xt3d, this%ja_xt3d, ierror)
169 call sparse_xt3d%destroy()
173 do i = 1, this%dis%nodes
175 do kk = this%ia_xt3d(i), this%ia_xt3d(i + 1) - 1
178 call sparse%addconnection(iglo, kglo, 1, iadded)
179 this%numextnbrs = this%numextnbrs + 1
184 call mem_allocate(this%ia_xt3d, 0,
'IA_XT3D', trim(this%memoryPath))
185 call mem_allocate(this%ja_xt3d, 0,
'JA_XT3D', trim(this%memoryPath))
196 integer(I4B),
intent(in) :: moffset
199 integer(I4B) :: i, j, iglo, jglo, niax, njax, ipos
200 integer(I4B) :: ipos_sln, icol_first, icol_last
201 integer(I4B) :: jj_xt3d
202 integer(I4B) :: igfirstnod, iglastnod
207 if (this%ixt3d == 1)
then
211 igfirstnod = moffset + 1
212 iglastnod = moffset + this%dis%nodes
215 niax = this%dis%nodes + 1
216 njax = this%numextnbrs
217 call mem_allocate(this%iax, niax,
'IAX', trim(this%memoryPath))
218 call mem_allocate(this%jax, njax,
'JAX', trim(this%memoryPath))
219 call mem_allocate(this%idxglox, njax,
'IDXGLOX', trim(this%memoryPath))
226 do i = 1, this%dis%nodes
230 icol_first = matrix_sln%get_first_col_pos(iglo)
231 icol_last = matrix_sln%get_last_col_pos(iglo)
232 do ipos_sln = icol_first, icol_last
236 jglo = matrix_sln%get_column(ipos_sln)
237 if (jglo < igfirstnod .or. jglo > iglastnod)
then
246 do jj_xt3d = this%ia_xt3d(i), this%ia_xt3d(i + 1) - 1
247 if (j == this%ja_xt3d(jj_xt3d))
then
255 this%jax(ipos) = matrix_sln%get_column(ipos_sln) - moffset
256 this%idxglox(ipos) = ipos_sln
261 this%iax(i + 1) = ipos
266 call mem_allocate(this%iax, 0,
'IAX', trim(this%memoryPath))
267 call mem_allocate(this%jax, 0,
'JAX', trim(this%memoryPath))
268 call mem_allocate(this%idxglox, 0,
'IDXGLOX', trim(this%memoryPath))
275 subroutine xt3d_ar(this, ibound, k11, ik33, k33, sat, ik22, k22, iangle1, &
276 iangle2, iangle3, angle1, angle2, angle3, inewton, icelltype)
281 integer(I4B),
dimension(:),
pointer,
contiguous,
intent(inout) :: ibound
282 real(DP),
dimension(:),
intent(in),
pointer,
contiguous :: k11
283 integer(I4B),
intent(in),
pointer :: ik33
284 real(DP),
dimension(:),
intent(in),
pointer,
contiguous :: k33
285 real(DP),
dimension(:),
intent(in),
pointer,
contiguous :: sat
286 integer(I4B),
intent(in),
pointer :: ik22
287 real(DP),
dimension(:),
intent(in),
pointer,
contiguous :: k22
288 integer(I4B),
intent(in),
pointer :: iangle1
289 integer(I4B),
intent(in),
pointer :: iangle2
290 integer(I4B),
intent(in),
pointer :: iangle3
291 real(DP),
dimension(:),
intent(in),
pointer,
contiguous :: angle1
292 real(DP),
dimension(:),
intent(in),
pointer,
contiguous :: angle2
293 real(DP),
dimension(:),
intent(in),
pointer,
contiguous :: angle3
294 integer(I4B),
intent(in),
pointer,
optional :: inewton
295 integer(I4B),
dimension(:),
intent(in),
pointer, &
296 contiguous,
optional :: icelltype
298 integer(I4B) :: n, nnbrs
300 character(len=*),
parameter :: fmtheader = &
301 "(1x, /1x, 'XT3D is active.'//)"
304 if (this%iout > 0)
then
305 write (this%iout, fmtheader)
309 this%ibound => ibound
316 this%iangle1 => iangle1
317 this%iangle2 => iangle2
318 this%iangle3 => iangle3
319 this%angle1 => angle1
320 this%angle2 => angle2
321 this%angle3 => angle3
323 if (
present(inewton))
then
325 this%inewton = inewton
327 if (
present(icelltype))
then
331 this%icelltype => icelltype
336 if (this%iangle2 == 0) this%nozee = .true.
340 do n = 1, this%dis%nodes
341 nnbrs = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
342 this%nbrmax = max(nnbrs, this%nbrmax)
346 if (this%dis%icondir == 0)
then
347 call store_error(
'Vertices not specified for discretization package, '// &
348 'but XT3D is active: '//trim(adjustl(this%memoryPath))// &
349 '. Vertices must be specified in discretization '// &
350 'package in order to use XT3D.', terminate=.true.)
354 if (this%dis%con%ianglex == 0)
then
355 call store_error(
'ANGLDEGX is not specified in the DIS package, '// &
356 'but XT3D is active: '//trim(adjustl(this%memoryPath))// &
357 '. ANGLDEGX must be provided in discretization '// &
358 'package in order to use XT3D.', terminate=.true.)
362 call this%allocate_arrays()
366 if (this%lamatsaved .and. .not. this%ldispersion) &
367 call this%xt3d_fcpc(this%dis%nodes, .true.)
371 subroutine xt3d_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
377 integer(I4B) :: kiter
379 integer(I4B),
intent(in),
dimension(:) :: idxglo
380 real(DP),
intent(inout),
dimension(:) :: rhs
381 real(DP),
intent(inout),
dimension(:) :: hnew
383 integer(I4B) :: nodes, nja
384 integer(I4B) :: n, m, ipos
385 logical :: allhc0, allhc1
386 integer(I4B) :: nnbr0, nnbr1
387 integer(I4B) :: il0, ii01, jjs01, il01, il10, ii00, ii11, ii10
389 integer(I4B),
dimension(this%nbrmax) :: inbr0, inbr1
390 real(DP) :: ar01, ar10
391 real(DP),
dimension(this%nbrmax, 3) :: vc0, vn0, vc1, vn1
392 real(DP),
dimension(this%nbrmax) :: dl0, dl0n, dl1, dl1n
393 real(DP),
dimension(3, 3) :: ck0, ck1
395 real(DP),
dimension(this%nbrmax) :: chati0, chat1j
396 real(DP) :: qnm, qnbrs
401 nodes = this%dis%nodes
402 nja = this%dis%con%nja
403 if (this%lamatsaved)
then
404 do i = 1, this%dis%con%nja
405 call matrix_sln%add_value_pos(idxglo(i), this%amatpc(i))
407 do i = 1, this%numextnbrs
408 call matrix_sln%add_value_pos(this%idxglox(i), this%amatpcx(i))
414 if (this%ibound(n) .eq. 0) cycle
416 if (this%lamatsaved)
then
417 if (this%iallpc(n) == 1) cycle
419 nnbr0 = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
421 call this%xt3d_load(nodes, n, nnbr0, inbr0, vc0, vn0, dl0, dl0n, &
426 ipos = this%dis%con%ia(n) + il0
427 if (this%dis%con%mask(ipos) == 0) cycle
431 if ((m .eq. 0) .or. (m .lt. n)) cycle
432 nnbr1 = this%dis%con%ia(m + 1) - this%dis%con%ia(m) - 1
434 call this%xt3d_load(nodes, m, nnbr1, inbr1, vc1, vn1, dl1, dl1n, &
437 call this%xt3d_indices(n, m, il0, ii01, jjs01, il01, il10, &
440 if (this%inewton /= 0)
then
444 call this%xt3d_areas(nodes, n, m, jjs01, .false., ar01, ar10, hnew)
448 call qconds(this%nbrmax, nnbr0, inbr0, il01, vc0, vn0, dl0, dl0n, ck0, &
449 nnbr1, inbr1, il10, vc1, vn1, dl1, dl1n, ck1, ar01, ar10, &
450 this%vcthresh, allhc0, allhc1, chat01, chati0, chat1j)
454 if (this%inewton /= 0)
then
456 qnm = chat01 * (hnew(m) - hnew(n))
458 call this%xt3d_qnbrs(nodes, n, m, nnbr0, inbr0, chati0, hnew, qnbrs)
461 call this%xt3d_qnbrs(nodes, m, n, nnbr1, inbr1, chat1j, hnew, qnbrs)
464 call this%xt3d_areas(nodes, n, m, jjs01, .true., ar01, ar10, hnew)
465 this%qsat(ii01) = qnm * ar01
467 call this%xt3d_areas(nodes, n, m, jjs01, .false., ar01, ar10, hnew)
468 chat01 = chat01 * ar01
469 chati0 = chati0 * ar01
470 chat1j = chat1j * ar01
474 call matrix_sln%add_value_pos(idxglo(ii00), -chat01)
475 call matrix_sln%add_value_pos(idxglo(ii01), chat01)
476 call matrix_sln%add_value_pos(idxglo(ii11), -chat01)
477 call matrix_sln%add_value_pos(idxglo(ii10), chat01)
478 if (this%ixt3d == 1)
then
479 call this%xt3d_amat_nbrs(nodes, n, ii00, nnbr0, nja, matrix_sln, &
480 inbr0, idxglo, chati0)
481 call this%xt3d_amat_nbrnbrs(nodes, n, m, ii01, nnbr1, nja, matrix_sln, &
482 inbr1, idxglo, chat1j)
483 call this%xt3d_amat_nbrs(nodes, m, ii11, nnbr1, nja, matrix_sln, &
484 inbr1, idxglo, chat1j)
485 call this%xt3d_amat_nbrnbrs(nodes, m, n, ii10, nnbr0, nja, matrix_sln, &
486 inbr0, idxglo, chati0)
488 call this%xt3d_rhs(nodes, n, m, nnbr0, inbr0, chati0, hnew, rhs)
489 call this%xt3d_rhs(nodes, m, n, nnbr1, inbr1, chat1j, hnew, rhs)
504 integer(I4B),
intent(in) :: nodes
505 logical,
intent(in) :: lsat
507 integer(I4B) :: n, m, ipos
509 logical :: allhc0, allhc1
510 integer(I4B) :: nnbr0, nnbr1
511 integer(I4B) :: il0, ii01, jjs01, il01, il10, ii00, ii11, ii10
512 integer(I4B),
dimension(this%nbrmax) :: inbr0, inbr1
513 real(DP) :: ar01, ar10
514 real(DP),
dimension(this%nbrmax, 3) :: vc0, vn0, vc1, vn1
515 real(DP),
dimension(this%nbrmax) :: dl0, dl0n, dl1, dl1n
516 real(DP),
dimension(3, 3) :: ck0, ck1
518 real(DP),
dimension(this%nbrmax) :: chati0, chat1j
521 do n = 1,
size(this%amatpc)
522 this%amatpc(n) =
dzero
524 do n = 1,
size(this%amatpcx)
525 this%amatpcx(n) =
dzero
532 if (this%iallpc(n) == 0) cycle
533 if (this%ibound(n) == 0) cycle
534 nnbr0 = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
536 call this%xt3d_load(nodes, n, nnbr0, inbr0, vc0, vn0, dl0, dl0n, &
541 ipos = this%dis%con%ia(n) + il0
542 if (this%dis%con%mask(ipos) == 0) cycle
547 nnbr1 = this%dis%con%ia(m + 1) - this%dis%con%ia(m) - 1
549 call this%xt3d_load(nodes, m, nnbr1, inbr1, vc1, vn1, dl1, dl1n, &
552 call this%xt3d_indices(n, m, il0, ii01, jjs01, il01, il10, &
555 call this%xt3d_areas(nodes, n, m, jjs01, lsat, ar01, ar10)
558 call qconds(this%nbrmax, nnbr0, inbr0, il01, vc0, vn0, dl0, dl0n, ck0, &
559 nnbr1, inbr1, il10, vc1, vn1, dl1, dl1n, ck1, ar01, ar10, &
560 this%vcthresh, allhc0, allhc1, chat01, chati0, chat1j)
562 this%amatpc(ii00) = this%amatpc(ii00) - chat01
563 this%amatpc(ii01) = this%amatpc(ii01) + chat01
564 this%amatpc(ii11) = this%amatpc(ii11) - chat01
565 this%amatpc(ii10) = this%amatpc(ii10) + chat01
566 call this%xt3d_amatpc_nbrs(nodes, n, ii00, nnbr0, inbr0, chati0)
567 call this%xt3d_amatpcx_nbrnbrs(nodes, n, m, ii01, nnbr1, inbr1, chat1j)
568 call this%xt3d_amatpc_nbrs(nodes, m, ii11, nnbr1, inbr1, chat1j)
569 call this%xt3d_amatpcx_nbrnbrs(nodes, m, n, ii10, nnbr0, inbr0, chati0)
576 subroutine xt3d_fhfb(this, kiter, nodes, nja, matrix_sln, idxglo, rhs, hnew, &
583 integer(I4B) :: kiter
584 integer(I4B),
intent(in) :: nodes
585 integer(I4B),
intent(in) :: nja
588 integer(I4B),
intent(in),
dimension(nja) :: idxglo
589 real(DP),
intent(inout),
dimension(nodes) :: rhs
590 real(DP),
intent(inout),
dimension(nodes) :: hnew
593 logical :: allhc0, allhc1
594 integer(I4B) :: nnbr0, nnbr1
595 integer(I4B) :: il0, ii01, jjs01, il01, il10, ii00, ii11, ii10, il
596 integer(I4B),
dimension(this%nbrmax) :: inbr0, inbr1
597 real(DP) :: ar01, ar10
598 real(DP),
dimension(this%nbrmax, 3) :: vc0, vn0, vc1, vn1
599 real(DP),
dimension(this%nbrmax) :: dl0, dl0n, dl1, dl1n
600 real(DP),
dimension(3, 3) :: ck0, ck1
602 real(DP),
dimension(this%nbrmax) :: chati0, chat1j
603 real(DP) :: qnm, qnbrs
609 nnbr0 = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
611 call this%xt3d_load(nodes, n, nnbr0, inbr0, vc0, vn0, dl0, dl0n, &
615 if (inbr0(il) .eq. m)
then
620 nnbr1 = this%dis%con%ia(m + 1) - this%dis%con%ia(m) - 1
622 call this%xt3d_load(nodes, m, nnbr1, inbr1, vc1, vn1, dl1, dl1n, &
625 call this%xt3d_indices(n, m, il0, ii01, jjs01, il01, il10, &
628 if (this%inewton /= 0)
then
632 call this%xt3d_areas(nodes, n, m, jjs01, .false., ar01, ar10, hnew)
637 call qconds(this%nbrmax, nnbr0, inbr0, il01, vc0, vn0, dl0, dl0n, &
638 ck0, nnbr1, inbr1, il10, vc1, vn1, dl1, dl1n, ck1, ar01, ar10, &
639 this%vcthresh, allhc0, allhc1, chat01, chati0, chat1j)
642 if (condhfb > dzero)
then
643 term = chat01 / (chat01 + condhfb)
647 chat01 = -chat01 * term
648 chati0 = -chati0 * term
649 chat1j = -chat1j * term
653 if (this%inewton /= 0)
then
655 qnm = chat01 * (hnew(m) - hnew(n))
657 call this%xt3d_qnbrs(nodes, n, m, nnbr0, inbr0, chati0, hnew, qnbrs)
660 call this%xt3d_qnbrs(nodes, m, n, nnbr1, inbr1, chat1j, hnew, qnbrs)
663 call this%xt3d_areas(nodes, n, m, jjs01, .true., ar01, ar10, hnew)
664 this%qsat(ii01) = this%qsat(ii01) + qnm * ar01
666 call this%xt3d_areas(nodes, n, m, jjs01, .false., ar01, ar10, hnew)
667 chat01 = chat01 * ar01
668 chati0 = chati0 * ar01
669 chat1j = chat1j * ar01
673 call matrix_sln%add_value_pos(idxglo(ii00), -chat01)
674 call matrix_sln%add_value_pos(idxglo(ii01), chat01)
675 call matrix_sln%add_value_pos(idxglo(ii11), -chat01)
676 call matrix_sln%add_value_pos(idxglo(ii10), chat01)
677 if (this%ixt3d == 1)
then
678 call this%xt3d_amat_nbrs(nodes, n, ii00, nnbr0, nja, matrix_sln, &
679 inbr0, idxglo, chati0)
680 call this%xt3d_amat_nbrnbrs(nodes, n, m, ii01, nnbr1, nja, matrix_sln, &
681 inbr1, idxglo, chat1j)
682 call this%xt3d_amat_nbrs(nodes, m, ii11, nnbr1, nja, matrix_sln, &
683 inbr1, idxglo, chat1j)
684 call this%xt3d_amat_nbrnbrs(nodes, m, n, ii10, nnbr0, nja, matrix_sln, &
685 inbr0, idxglo, chati0)
687 call this%xt3d_rhs(nodes, n, m, nnbr0, inbr0, chati0, hnew, rhs)
688 call this%xt3d_rhs(nodes, m, n, nnbr1, inbr1, chat1j, hnew, rhs)
693 subroutine xt3d_fn(this, kiter, nodes, nja, matrix_sln, idxglo, rhs, hnew)
699 integer(I4B) :: kiter
700 integer(I4B),
intent(in) :: nodes
701 integer(I4B),
intent(in) :: nja
703 integer(I4B),
intent(in),
dimension(nja) :: idxglo
704 real(DP),
intent(inout),
dimension(nodes) :: rhs
705 real(DP),
intent(inout),
dimension(nodes) :: hnew
707 integer(I4B) :: n, m, ipos
708 integer(I4B) :: nnbr0
709 integer(I4B) :: il0, ii01, jjs01, il01, il10, ii00, ii11, ii10
710 integer(I4B),
dimension(this%nbrmax) :: inbr0
711 integer(I4B) :: iups, idn
712 real(DP) :: topup, botup, derv, term
718 if (this%ibound(n) .eq. 0) cycle
722 if (this%lamatsaved)
then
723 if (this%iallpc(n) == 1) cycle
725 nnbr0 = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
729 call this%xt3d_load_inbr(n, nnbr0, inbr0)
734 ipos = this%dis%con%ia(n) + il0
735 if (this%dis%con%mask(ipos) == 0) cycle
740 if ((inbr0(il0) .eq. 0) .or. (m .lt. n)) cycle
743 call this%xt3d_indices(n, m, il0, ii01, jjs01, il01, il10, &
748 if (hnew(m) < hnew(n)) iups = n
750 if (iups == n) idn = m
753 if ((this%icelltype(iups) == 0) .and. (this%ixt3d .eq. 1)) cycle
757 topup = this%dis%top(iups)
758 botup = this%dis%bot(iups)
759 if (this%dis%con%ihc(jjs01) == 2)
then
760 topup = min(this%dis%top(n), this%dis%top(m))
761 botup = max(this%dis%bot(n), this%dis%bot(m))
766 term = this%qsat(ii01) * derv
772 call matrix_sln%add_value_pos(idxglo(ii00), term)
773 rhs(n) = rhs(n) + term * hnew(n)
776 call matrix_sln%add_value_pos(idxglo(ii10), -term)
777 rhs(m) = rhs(m) - term * hnew(n)
783 call matrix_sln%add_value_pos(idxglo(ii01), term)
784 rhs(n) = rhs(n) + term * hnew(m)
787 call matrix_sln%add_value_pos(idxglo(ii11), -term)
788 rhs(m) = rhs(m) - term * hnew(m)
802 real(DP),
intent(inout),
dimension(:) :: hnew
803 real(DP),
intent(inout),
dimension(:) :: flowja
805 integer(I4B) :: n, ipos, m, nodes
806 real(DP) :: qnm, qnbrs
807 logical :: allhc0, allhc1
808 integer(I4B) :: nnbr0, nnbr1
809 integer(I4B) :: il0, ii01, jjs01, il01, il10, ii00, ii11, ii10
810 integer(I4B),
dimension(this%nbrmax) :: inbr0, inbr1
811 real(DP) :: ar01, ar10
812 real(DP),
dimension(this%nbrmax, 3) :: vc0, vn0, vc1, vn1
813 real(DP),
dimension(this%nbrmax) :: dl0, dl0n, dl1, dl1n
814 real(DP),
dimension(3, 3) :: ck0, ck1
816 real(DP),
dimension(this%nbrmax) :: chati0, chat1j
819 nodes = this%dis%nodes
823 if (this%ibound(n) .eq. 0) cycle
824 nnbr0 = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
827 call this%xt3d_load(nodes, n, nnbr0, inbr0, vc0, vn0, dl0, dl0n, &
836 if ((inbr0(il0) .eq. 0) .or. (m .lt. n)) cycle
837 nnbr1 = this%dis%con%ia(m + 1) - this%dis%con%ia(m) - 1
840 call this%xt3d_load(nodes, m, nnbr1, inbr1, vc1, vn1, dl1, dl1n, &
844 call this%xt3d_indices(n, m, il0, ii01, jjs01, il01, il10, &
848 if (this%inewton /= 0) &
849 call this%xt3d_areas(nodes, n, m, jjs01, .true., ar01, ar10, hnew)
850 call this%xt3d_areas(nodes, n, m, jjs01, .false., ar01, ar10, hnew)
854 call qconds(this%nbrmax, nnbr0, inbr0, il01, vc0, vn0, dl0, dl0n, ck0, &
855 nnbr1, inbr1, il10, vc1, vn1, dl1, dl1n, ck1, ar01, ar10, &
856 this%vcthresh, allhc0, allhc1, chat01, chati0, chat1j)
859 qnm = chat01 * (hnew(m) - hnew(n))
862 call this%xt3d_qnbrs(nodes, n, m, nnbr0, inbr0, chati0, hnew, qnbrs)
866 call this%xt3d_qnbrs(nodes, m, n, nnbr1, inbr1, chat1j, hnew, qnbrs)
869 flowja(ipos) = flowja(ipos) + qnm
870 flowja(this%dis%con%isym(ipos)) = flowja(this%dis%con%isym(ipos)) - qnm
884 real(DP),
intent(inout),
dimension(:) :: hnew
885 real(DP),
intent(inout),
dimension(:) :: flowja
888 integer(I4B) :: nodes
889 logical :: allhc0, allhc1
891 integer(I4B) :: nnbr0, nnbr1
892 integer(I4B) :: il0, ii01, jjs01, il01, il10, ii00, ii11, ii10, il
893 integer(I4B),
dimension(this%nbrmax) :: inbr0, inbr1
895 real(DP) :: ar01, ar10
896 real(DP),
dimension(this%nbrmax, 3) :: vc0, vn0, vc1, vn1
897 real(DP),
dimension(this%nbrmax) :: dl0, dl0n, dl1, dl1n
898 real(DP),
dimension(3, 3) :: ck0, ck1
900 real(DP),
dimension(this%nbrmax) :: chati0, chat1j
901 real(DP) :: qnm, qnbrs
906 nodes = this%dis%nodes
907 nnbr0 = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
910 call this%xt3d_load(nodes, n, nnbr0, inbr0, vc0, vn0, dl0, dl0n, &
915 if (inbr0(il) .eq. m)
then
921 nnbr1 = this%dis%con%ia(m + 1) - this%dis%con%ia(m) - 1
924 call this%xt3d_load(nodes, m, nnbr1, inbr1, vc1, vn1, dl1, dl1n, &
928 call this%xt3d_indices(n, m, il0, ii01, jjs01, il01, il10, &
932 if (this%inewton /= 0)
then
936 call this%xt3d_areas(nodes, n, m, jjs01, .false., ar01, ar10, hnew)
941 call qconds(this%nbrmax, nnbr0, inbr0, il01, vc0, vn0, dl0, dl0n, &
942 ck0, nnbr1, inbr1, il10, vc1, vn1, dl1, dl1n, ck1, ar01, ar10, &
943 this%vcthresh, allhc0, allhc1, chat01, chati0, chat1j)
946 if (condhfb > dzero)
then
947 term = chat01 / (chat01 + condhfb)
951 chat01 = -chat01 * term
952 chati0 = -chati0 * term
953 chat1j = -chat1j * term
956 qnm = chat01 * (hnew(m) - hnew(n))
959 call this%xt3d_qnbrs(nodes, n, m, nnbr0, inbr0, chati0, hnew, qnbrs)
963 call this%xt3d_qnbrs(nodes, m, n, nnbr1, inbr1, chat1j, hnew, qnbrs)
968 if (this%inewton /= 0)
then
969 call this%xt3d_areas(nodes, n, m, jjs01, .true., ar01, ar10, hnew)
970 call this%xt3d_areas(nodes, n, m, jjs01, .false., ar01, ar10, hnew)
975 flowja(ipos) = flowja(ipos) + qnm
976 flowja(this%dis%con%isym(ipos)) = flowja(this%dis%con%isym(ipos)) - qnm
988 if (this%ixt3d /= 0)
then
1023 call mem_allocate(this%ixt3d,
'IXT3D', this%memoryPath)
1024 call mem_allocate(this%nbrmax,
'NBRMAX', this%memoryPath)
1025 call mem_allocate(this%inunit,
'INUNIT', this%memoryPath)
1027 call mem_allocate(this%inewton,
'INEWTON', this%memoryPath)
1028 call mem_allocate(this%numextnbrs,
'NUMEXTNBRS', this%memoryPath)
1029 call mem_allocate(this%nozee,
'NOZEE', this%memoryPath)
1030 call mem_allocate(this%vcthresh,
'VCTHRESH', this%memoryPath)
1031 call mem_allocate(this%lamatsaved,
'LAMATSAVED', this%memoryPath)
1032 call mem_allocate(this%ldispersion,
'LDISPERSION', this%memoryPath)
1041 this%nozee = .false.
1042 this%vcthresh = 1.d-10
1043 this%lamatsaved = .false.
1044 this%ldispersion = .false.
1055 integer(I4B) :: njax
1059 if (this%inewton /= 0)
then
1060 call mem_allocate(this%qsat, this%dis%nja,
'QSAT', this%memoryPath)
1062 call mem_allocate(this%qsat, 0,
'QSAT', this%memoryPath)
1068 if (this%ldispersion)
then
1072 this%lamatsaved = .true.
1073 call mem_allocate(this%iallpc, this%dis%nodes,
'IALLPC', this%memoryPath)
1074 do n = 1, this%dis%nodes
1082 call this%xt3d_iallpc()
1086 if (this%lamatsaved)
then
1087 call mem_allocate(this%amatpc, this%dis%nja,
'AMATPC', this%memoryPath)
1088 njax = this%numextnbrs
1089 call mem_allocate(this%amatpcx, njax,
'AMATPCX', this%memoryPath)
1091 call mem_allocate(this%amatpc, 0,
'AMATPC', this%memoryPath)
1092 call mem_allocate(this%amatpcx, 0,
'AMATPCX', this%memoryPath)
1096 call mem_allocate(this%rmatck, 3, 3,
'RMATCK', this%memoryPath)
1100 if (this%inewton /= 0)
then
1102 else if (this%lamatsaved)
then
1104 this%amatpcx =
dzero
1115 integer(I4B) :: n, m, mm, il0, il1
1116 integer(I4B) :: nnbr0, nnbr1
1117 integer(I4B),
dimension(this%nbrmax) :: inbr0, inbr1
1119 if (this%ixt3d == 2)
then
1120 this%lamatsaved = .false.
1121 call mem_allocate(this%iallpc, 0,
'IALLPC', this%memoryPath)
1125 call mem_allocate(this%iallpc, this%dis%nodes,
'IALLPC', this%memoryPath)
1126 do n = 1, this%dis%nodes
1132 do n = 1, this%dis%nodes
1133 if (this%icelltype(n) /= 0)
then
1137 nnbr0 = this%dis%con%ia(n + 1) - this%dis%con%ia(n) - 1
1138 call this%xt3d_load_inbr(n, nnbr0, inbr0)
1142 if (this%icelltype(m) /= 0)
then
1147 nnbr1 = this%dis%con%ia(m + 1) - this%dis%con%ia(m) - 1
1148 call this%xt3d_load_inbr(m, nnbr1, inbr1)
1151 if (this%icelltype(mm) /= 0)
then
1163 this%lamatsaved = .false.
1164 do n = 1, this%dis%nodes
1165 if (this%iallpc(n) == 1)
then
1166 this%lamatsaved = .true.
1172 if (.not. this%lamatsaved)
then
1185 integer(I4B) :: n, m, il0, ii01, jjs01, il01, il10, ii00, ii11, ii10
1187 integer(I4B) :: iinm
1194 call this%xt3d_get_iinm(m, n, iinm)
1195 il10 = iinm - this%dis%con%ia(m)
1197 ii00 = this%dis%con%ia(n)
1201 jjs01 = this%dis%con%jas(ii01)
1203 ii11 = this%dis%con%ia(m)
1211 subroutine xt3d_load(this, nodes, n, nnbr, inbr, vc, vn, dl, dln, ck, allhc)
1217 integer(I4B),
intent(in) :: nodes
1218 integer(I4B) :: n, nnbr
1219 integer(I4B),
dimension(this%nbrmax) :: inbr
1220 real(DP),
dimension(this%nbrmax, 3) :: vc, vn
1221 real(DP),
dimension(this%nbrmax) :: dl, dln
1222 real(DP),
dimension(3, 3) :: ck
1224 integer(I4B) :: il, ii, jj, jjs
1225 integer(I4B) :: ihcnjj
1226 real(DP) :: satn, satjj
1227 real(DP) :: cl1njj, cl2njj, dltot, ooclsum
1231 ck(1, 1) = this%k11(n)
1232 ck(2, 2) = this%k22(n)
1233 ck(3, 3) = this%k33(n)
1234 call this%xt3d_fillrmatck(n)
1235 ck = matmul(this%rmatck, ck)
1236 ck = matmul(ck, transpose(this%rmatck))
1244 ii = il + this%dis%con%ia(n)
1245 jj = this%dis%con%ja(ii)
1246 jjs = this%dis%con%jas(ii)
1247 if (this%ibound(jj) .ne. 0)
then
1250 satjj = this%sat(jj)
1253 ihcnjj = this%dis%con%ihc(jjs)
1254 call this%dis%connection_normal(n, jj, ihcnjj, vn(il, 1), vn(il, 2), &
1256 call this%dis%connection_vector(n, jj, this%nozee, satn, satjj, ihcnjj, &
1257 vc(il, 1), vc(il, 2), vc(il, 3), dltot)
1259 cl1njj = this%dis%con%cl1(jjs)
1260 cl2njj = this%dis%con%cl2(jjs)
1262 cl1njj = this%dis%con%cl2(jjs)
1263 cl2njj = this%dis%con%cl1(jjs)
1265 ooclsum = 1d0 / (cl1njj + cl2njj)
1266 dl(il) = dltot * cl1njj * ooclsum
1267 dln(il) = dltot * cl2njj * ooclsum
1268 if (this%dis%con%ihc(jjs) .eq. 0) allhc = .false.
1280 integer(I4B) :: n, nnbr
1281 integer(I4B),
dimension(this%nbrmax) :: inbr
1283 integer(I4B) :: il, ii, jj
1288 ii = il + this%dis%con%ia(n)
1289 jj = this%dis%con%ja(ii)
1290 if (this%ibound(jj) .ne. 0)
then
1300 subroutine xt3d_areas(this, nodes, n, m, jjs01, lsat, ar01, ar10, hnew)
1304 integer(I4B) :: nodes, n, m, jjs01
1305 real(DP) :: ar01, ar10
1306 real(DP),
intent(inout),
dimension(:),
optional :: hnew
1308 real(DP) :: topn, botn, topm, botm, thksatn, thksatm
1309 real(DP) :: sill_top, sill_bot, tpn, tpm
1313 if (this%dis%con%ihc(jjs01) == 0)
then
1316 ar01 = this%dis%con%hwva(jjs01)
1318 else if (this%inewton /= 0)
then
1324 topn = this%dis%top(n)
1325 botn = this%dis%bot(n)
1326 topm = this%dis%top(m)
1327 botm = this%dis%bot(m)
1328 thksatn = topn - botn
1329 thksatm = topm - botm
1330 if (this%dis%con%ihc(jjs01) .eq. 2)
then
1332 sill_top = min(topn, topm)
1333 sill_bot = max(botn, botm)
1334 tpn = botn + thksatn
1335 tpm = botm + thksatm
1336 thksatn = max(min(tpn, sill_top) - sill_bot,
dzero)
1337 thksatm = max(min(tpm, sill_top) - sill_bot,
dzero)
1339 ar01 = this%dis%con%hwva(jjs01) *
dhalf * (thksatn + thksatm)
1345 if (hnew(m) < hnew(n))
then
1346 satups = this%sat(n)
1348 satups = this%sat(m)
1350 ar01 = ar01 * satups
1356 topn = this%dis%top(n)
1357 botn = this%dis%bot(n)
1358 topm = this%dis%top(m)
1359 botm = this%dis%bot(m)
1360 thksatn = topn - botn
1361 thksatm = topm - botm
1362 if (.not. lsat)
then
1363 thksatn = this%sat(n) * thksatn
1364 thksatm = this%sat(m) * thksatm
1366 if (this%dis%con%ihc(jjs01) == 2)
then
1368 sill_top = min(topn, topm)
1369 sill_bot = max(botn, botm)
1370 tpn = botn + thksatn
1371 tpm = botm + thksatm
1372 thksatn = max(min(tpn, sill_top) - sill_bot,
dzero)
1373 thksatm = max(min(tpm, sill_top) - sill_bot,
dzero)
1375 ar01 = this%dis%con%hwva(jjs01) * thksatn
1376 ar10 = this%dis%con%hwva(jjs01) * thksatm
1383 matrix_sln, inbr, idxglo, chat)
1386 integer(I4B),
intent(in) :: nodes
1387 integer(I4B) :: n, idiag, nnbr, nja
1389 integer(I4B),
dimension(this%nbrmax) :: inbr
1390 integer(I4B),
intent(in),
dimension(nja) :: idxglo
1391 real(DP),
dimension(this%nbrmax) :: chat
1393 integer(I4B) :: iil, iii
1396 if (inbr(iil) .ne. 0)
then
1397 iii = this%dis%con%ia(n) + iil
1398 call matrix_sln%add_value_pos(idxglo(idiag), -chat(iil))
1399 call matrix_sln%add_value_pos(idxglo(iii), chat(iil))
1407 matrix_sln, inbr, idxglo, chat)
1410 integer(I4B),
intent(in) :: nodes
1411 integer(I4B) :: n, m, ii01, nnbr, nja
1413 integer(I4B),
dimension(this%nbrmax) :: inbr
1414 integer(I4B),
intent(in),
dimension(nja) :: idxglo
1415 real(DP),
dimension(this%nbrmax) :: chat
1417 integer(I4B) :: iil, iii, jjj, iixjjj, iijjj
1420 if (inbr(iil) .ne. 0)
then
1421 call matrix_sln%add_value_pos(idxglo(ii01), chat(iil))
1422 iii = this%dis%con%ia(m) + iil
1423 jjj = this%dis%con%ja(iii)
1424 call this%xt3d_get_iinmx(n, jjj, iixjjj)
1425 if (iixjjj .ne. 0)
then
1426 call matrix_sln%add_value_pos(this%idxglox(iixjjj), -chat(iil))
1428 call this%xt3d_get_iinm(n, jjj, iijjj)
1429 call matrix_sln%add_value_pos(idxglo(iijjj), -chat(iil))
1440 integer(I4B),
intent(in) :: nodes
1441 integer(I4B) :: n, idiag, nnbr
1442 integer(I4B),
dimension(this%nbrmax) :: inbr
1443 real(DP),
dimension(this%nbrmax) :: chat
1445 integer(I4B) :: iil, iii
1448 iii = this%dis%con%ia(n) + iil
1449 this%amatpc(idiag) = this%amatpc(idiag) - chat(iil)
1450 this%amatpc(iii) = this%amatpc(iii) + chat(iil)
1459 integer(I4B),
intent(in) :: nodes
1460 integer(I4B) :: n, m, ii01, nnbr
1461 integer(I4B),
dimension(this%nbrmax) :: inbr
1462 real(DP),
dimension(this%nbrmax) :: chat
1464 integer(I4B) :: iil, iii, jjj, iixjjj, iijjj
1467 this%amatpc(ii01) = this%amatpc(ii01) + chat(iil)
1468 iii = this%dis%con%ia(m) + iil
1469 jjj = this%dis%con%ja(iii)
1470 call this%xt3d_get_iinmx(n, jjj, iixjjj)
1471 if (iixjjj .ne. 0)
then
1472 this%amatpcx(iixjjj) = this%amatpcx(iixjjj) - chat(iil)
1474 call this%xt3d_get_iinm(n, jjj, iijjj)
1475 this%amatpc(iijjj) = this%amatpc(iijjj) - chat(iil)
1486 integer(I4B) :: n, m, iinm
1488 integer(I4B) :: ii, jj
1491 do ii = this%dis%con%ia(n), this%dis%con%ia(n + 1) - 1
1492 jj = this%dis%con%ja(ii)
1506 integer(I4B) :: n, m, iinmx
1508 integer(I4B) :: iix, jjx
1511 do iix = this%iax(n), this%iax(n + 1) - 1
1513 if (jjx .eq. m)
then
1522 subroutine xt3d_rhs(this, nodes, n, m, nnbr, inbr, chat, hnew, &
1526 integer(I4B),
intent(in) :: nodes
1527 integer(I4B) :: n, m, nnbr
1528 integer(I4B),
dimension(this%nbrmax) :: inbr
1529 real(DP),
dimension(this%nbrmax) :: chat
1530 real(DP),
intent(inout),
dimension(nodes) :: hnew, rhs
1532 integer(I4B) :: iil, iii, jjj
1536 if (inbr(iil) .ne. 0)
then
1537 iii = iil + this%dis%con%ia(n)
1538 jjj = this%dis%con%ja(iii)
1539 term = chat(iil) * (hnew(jjj) - hnew(n))
1540 rhs(n) = rhs(n) - term
1541 rhs(m) = rhs(m) + term
1548 subroutine xt3d_qnbrs(this, nodes, n, m, nnbr, inbr, chat, hnew, &
1552 integer(I4B),
intent(in) :: nodes
1553 integer(I4B) :: n, m, nnbr
1554 integer(I4B),
dimension(this%nbrmax) :: inbr
1556 real(DP),
dimension(this%nbrmax) :: chat
1557 real(DP),
intent(inout),
dimension(nodes) :: hnew
1559 integer(I4B) :: iil, iii, jjj
1564 if (inbr(iil) .ne. 0)
then
1565 iii = iil + this%dis%con%ia(n)
1566 jjj = this%dis%con%ja(iii)
1567 term = chat(iil) * (hnew(jjj) - hnew(n))
1568 qnbrs = qnbrs + term
1580 integer(I4B),
intent(in) :: n
1582 real(DP) :: ang1, ang2, ang3, ang2d, ang3d
1583 real(DP) :: s1, c1, s2, c2, s3, c3
1585 if (this%nozee)
then
1588 ang1 = this%angle1(n)
1592 ang1 = this%angle1(n)
1593 ang2 = this%angle2(n)
1594 ang3 = this%angle3(n)
1602 this%rmatck(1, 1) = c1 * c2
1603 this%rmatck(1, 2) = c1 * s2 * s3 - s1 * c3
1604 this%rmatck(1, 3) = -c1 * s2 * c3 - s1 * s3
1605 this%rmatck(2, 1) = s1 * c2
1606 this%rmatck(2, 2) = s1 * s2 * s3 + c1 * c3
1607 this%rmatck(2, 3) = -s1 * s2 * c3 + c1 * s3
1608 this%rmatck(3, 1) = s2
1609 this%rmatck(3, 2) = -c2 * s3
1610 this%rmatck(3, 3) = c2 * c3
This module contains simulation constants.
real(dp), parameter dhalf
real constant 1/2
real(dp), parameter dzero
real constant zero
integer(i4b), parameter lenmempath
maximum length of the memory path
real(dp), parameter done
real constant 1
This module defines variable data types.
character(len=lenmempath) function create_mem_path(component, subcomponent, context)
returns the path to the memory object
This module contains simulation methods.
subroutine, public store_error(msg, terminate)
Store an error message.
real(dp) function squadraticsaturationderivative(top, bot, x, eps)
@ brief Derivative of the quadratic saturation function
subroutine qconds(nnbrmx, nnbr0, inbr0, il01, vc0, vn0, dl0, dl0n, ck0, nnbr1, inbr1, il10, vc1, vn1, dl1, dl1n, ck1, ar01, ar10, vcthresh, allhc0, allhc1, chat01, chati0, chat1j)
Compute the "conductances" in the normal-flux expression for an interface (modflow-usg version)....
subroutine xt3d_load(this, nodes, n, nnbr, inbr, vc, vn, dl, dln, ck, allhc)
Load conductivity and connection info for a cell into arrays used by XT3D.
subroutine xt3d_ac(this, moffset, sparse)
Add connections for extended neighbors to the sparse matrix.
subroutine xt3d_get_iinmx(this, n, m, iinmx)
Get position of n-m "extended connection" in jax array (return 0 if not connected)
subroutine xt3d_amatpcx_nbrnbrs(this, nodes, n, m, ii01, nnbr, inbr, chat)
Add contributions from neighbors of neighbor to amatpc and amatpcx.
subroutine allocate_arrays(this)
Allocate xt3d arrays.
subroutine xt3d_amat_nbrnbrs(this, nodes, n, m, ii01, nnbr, nja, matrix_sln, inbr, idxglo, chat)
Add contributions from neighbors of neighbor to amat.
subroutine xt3d_da(this)
Deallocate variables.
subroutine allocate_scalars(this)
Allocate scalar pointer variables.
subroutine xt3d_amat_nbrs(this, nodes, n, idiag, nnbr, nja, matrix_sln, inbr, idxglo, chat)
Add contributions from neighbors to amat.
subroutine xt3d_areas(this, nodes, n, m, jjs01, lsat, ar01, ar10, hnew)
Compute interfacial areas.
subroutine xt3d_flowja(this, hnew, flowja)
Budget.
subroutine xt3d_get_iinm(this, n, m, iinm)
Get position of n-m connection in ja array (return 0 if not connected)
subroutine xt3d_fhfb(this, kiter, nodes, nja, matrix_sln, idxglo, rhs, hnew, n, m, condhfb)
Formulate HFB correction.
subroutine xt3d_flowjahfb(this, n, m, hnew, flowja, condhfb)
hfb contribution to flowja when xt3d is used
subroutine xt3d_indices(this, n, m, il0, ii01, jjs01, il01, il10, ii00, ii11, ii10)
Set various indices for XT3D.
subroutine xt3d_load_inbr(this, n, nnbr, inbr)
Load neighbor list for a cell.
subroutine xt3d_fcpc(this, nodes, lsat)
Formulate for permanently confined connections and save in amatpc and amatpcx.
subroutine xt3d_fn(this, kiter, nodes, nja, matrix_sln, idxglo, rhs, hnew)
Fill Newton terms for xt3d.
subroutine xt3d_fillrmatck(this, n)
Fill rmat array for cell n.
subroutine xt3d_mc(this, moffset, matrix_sln)
Map connections and construct iax, jax, and idxglox.
subroutine, public xt3d_cr(xt3dobj, name_model, inunit, iout, ldispopt)
Create a new xt3d object.
subroutine xt3d_qnbrs(this, nodes, n, m, nnbr, inbr, chat, hnew, qnbrs)
Add contributions to flow from neighbors.
subroutine xt3d_fc(this, kiter, matrix_sln, idxglo, rhs, hnew)
Formulate.
subroutine xt3d_rhs(this, nodes, n, m, nnbr, inbr, chat, hnew, rhs)
Add contributions to rhs.
subroutine xt3d_df(this, dis)
Define the xt3d object.
subroutine xt3d_iallpc(this)
Allocate and populate iallpc array. Set lamatsaved.
subroutine xt3d_ar(this, ibound, k11, ik33, k33, sat, ik22, k22, iangle1, iangle2, iangle3, angle1, angle2, angle3, inewton, icelltype)
Allocate and Read.
subroutine xt3d_amatpc_nbrs(this, nodes, n, idiag, nnbr, inbr, chat)
Add contributions from neighbors to amatpc.