My Project
Public Member Functions | List of all members
mod_utils::interp_pnodal Interface Reference

Public Member Functions

real(spa) function interp_pnodal_2d_flt (xloc, yloc, i, Field)
 
real(spa) function interp_pnodal_3d_flt (xloc, yloc, sigloc, lvls, i, Field)
 
real(dp) function interp_pnodal_2d_dbl (xloc, yloc, i, Field)
 
real(dp) function interp_pnodal_3d_dbl (xloc, yloc, sigloc, lvls, i, Field)
 

Detailed Description

Definition at line 114 of file mod_utils.f90.

Member Function/Subroutine Documentation

◆ interp_pnodal_2d_dbl()

real(dp) function mod_utils::interp_pnodal::interp_pnodal_2d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
integer, intent(in)  i,
real(dp), dimension(:), intent(in), pointer  Field 
)

Definition at line 585 of file mod_utils.f90.

585  USE all_vars, only : aw0,awx,awy,nv,xc,yc
586  IMPLICIT NONE
587  REAL(DP) :: FPT
588  INTEGER, INTENT(IN) :: i ! Cell Number
589  REAL(DP), INTENT(IN):: xloc,yloc
590  REAL(DP), POINTER, INTENT(IN) :: FIELD(:)
591 
592  REAL(DP):: X0c, Y0c,F0,Fx,Fy
593  INTEGER :: n1,n2,n3
594  REAL(DP) :: dx_sph
595  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PNODAL_2D"
596 
597  !element location (i) and surrounding nodes (n1,n2,n3)
598  n1 = nv(i,1)
599  n2 = nv(i,2)
600  n3 = nv(i,3)
601 
602 !offset from element center
603 !---------------------------------------------------------------
604  !offset from element center
605  x0c = xloc - xc(i)
606  y0c = yloc - yc(i)
607 !---------------------------------------------------------------
608 !---------------------------------------------------------------
609 
610  !linear interpolation of Field
611  f0 = aw0(i,1)*field(n1)+aw0(i,2)*field(n2)+aw0(i,3)*field(n3)
612  fx = awx(i,1)*field(n1)+awx(i,2)*field(n2)+awx(i,3)*field(n3)
613  fy = awy(i,1)*field(n1)+awy(i,2)*field(n2)+awy(i,3)*field(n3)
614  fpt = f0 + fx*x0c + fy*y0c
615 
616 
617  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"END: INTERP_PNODAL_2D"
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target awx
Definition: mod_main.f90:1333
real(sp), dimension(:,:), allocatable, target aw0
Definition: mod_main.f90:1335
real(sp), dimension(:,:), allocatable, target awy
Definition: mod_main.f90:1334
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003

◆ interp_pnodal_2d_flt()

real(spa) function mod_utils::interp_pnodal::interp_pnodal_2d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
integer, intent(in)  i,
real(spa), dimension(:), intent(in), pointer  Field 
)

Definition at line 548 of file mod_utils.f90.

548  USE all_vars, only : aw0,awx,awy,nv,xc,yc
549  IMPLICIT NONE
550  REAL(SPA) :: FPT
551  INTEGER, INTENT(IN) :: i ! Cell Number
552  REAL(SPA), INTENT(IN):: xloc,yloc
553  REAL(SPA), POINTER, INTENT(IN) :: FIELD(:)
554 
555  REAL(SPA):: X0c, Y0c,F0,Fx,Fy
556  INTEGER :: n1,n2,n3
557  REAL(SPA) :: dx_sph
558  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PNODAL_2D"
559 
560  !element location (i) and surrounding nodes (n1,n2,n3)
561  n1 = nv(i,1)
562  n2 = nv(i,2)
563  n3 = nv(i,3)
564 
565 !offset from element center
566 !---------------------------------------------------------------
567 
568  x0c = xloc - xc(i)
569  y0c = yloc - yc(i)
570 !---------------------------------------------------------------
571 !---------------------------------------------------------------
572 
573  !linear interpolation of Field
574  f0 = aw0(i,1)*field(n1)+aw0(i,2)*field(n2)+aw0(i,3)*field(n3)
575  fx = awx(i,1)*field(n1)+awx(i,2)*field(n2)+awx(i,3)*field(n3)
576  fy = awy(i,1)*field(n1)+awy(i,2)*field(n2)+awy(i,3)*field(n3)
577  fpt = f0 + fx*x0c + fy*y0c
578 
579 
580  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"END: INTERP_PNODAL_2D"
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target awx
Definition: mod_main.f90:1333
real(sp), dimension(:,:), allocatable, target aw0
Definition: mod_main.f90:1335
real(sp), dimension(:,:), allocatable, target awy
Definition: mod_main.f90:1334
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003

◆ interp_pnodal_3d_dbl()

real(dp) function mod_utils::interp_pnodal::interp_pnodal_3d_dbl ( real(dp), intent(in)  xloc,
real(dp), intent(in)  yloc,
real(dp), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(dp), dimension(:,:), intent(in), pointer  Field 
)

Definition at line 752 of file mod_utils.f90.

752  USE all_vars, only : aw0,awx,awy,nv,xc,yc,kbm2,kbm1,kb,z1,zz1,dz1,dzz1
753  IMPLICIT NONE
754  REAL(DP) :: FPT
755  INTEGER, INTENT(IN) :: i, lvls ! Cell Number, Number of Levels in data (kb or kbm1)
756  REAL(DP), INTENT(IN):: xloc,yloc,sigloc
757  REAL(DP), POINTER, INTENT(IN) :: FIELD(:,:)
758 
759  REAL(DP):: X0c,Y0c,F0,Fx,Fy,F_LOWER,F_UPPER, alpha,dsig
760  INTEGER :: n1,n2,n3,k1,k2,k
761  REAL(DP) :: dx_sph
762 
763  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PNODAL_3D"
764 
765  !element location (i) and surrounding nodes (n1,n2,n3)
766  n1 = nv(i,1)
767  n2 = nv(i,2)
768  n3 = nv(i,3)
769 
770 !offset from element center
771 !---------------------------------------------------------------
772 
773  x0c = xloc - xc(i)
774  y0c = yloc - yc(i)
775 !---------------------------------------------------------------
776 !---------------------------------------------------------------
777 
778 
779  ! Determine the layer in which the point resides
780  IF(lvls == kbm1) THEN
781 
782  !top
783  if(sigloc >= zz1(i,1))then
784  k1 = 1
785  k2 = 1
786  alpha = -1
787  elseif(sigloc > zz1(i,kbm1)) then!intermediate
788  do k=1,kbm2
789  if(sigloc < zz1(i,k) .and. sigloc >= zz1(i,k+1) )then
790  k1 = k
791  k2 = k+1
792  alpha = (zz1(i,k)-sigloc)/dzz1(i,k)
793  exit
794  endif
795  end do
796  else
797  ! TOP
798  k1 = kbm1
799  k2 = kbm1
800  alpha = -1
801  endif
802  ELSE IF(lvls == kb) THEN
803 
804  !top
805  if(sigloc >= z1(i,1))then
806  k1 = 1
807  k2 = 1
808  alpha = -1
809  elseif(sigloc > z1(i,kb)) then !intermediate
810  do k=1,kbm1
811  if(sigloc < z1(i,k) .and. sigloc >= z1(i,k+1) )then
812  k1 = k
813  k2 = k+1
814  alpha = (z1(i,k)-sigloc)/dz1(i,k)
815  endif
816  end do
817  else
818  !bottom
819  k1 = kbm1
820  k2 = kbm1
821  alpha = -1
822  endif
823 
824  ELSE
825  CALL fatal_error("INTERP_PNODAL_3D: Invalid number of levels passed",&
826  & "(Must be equal to either KB or KBM1")
827  END IF
828 
829  !linear interpolation of Field
830  f0 = aw0(i,1)*field(n1,k1)+aw0(i,2)*field(n2,k1)+aw0(i,3)*field(n3,k1)
831  fx = awx(i,1)*field(n1,k1)+awx(i,2)*field(n2,k1)+awx(i,3)*field(n3,k1)
832  fy = awy(i,1)*field(n1,k1)+awy(i,2)*field(n2,k1)+awy(i,3)*field(n3,k1)
833  f_upper = f0 + fx*x0c + fy*y0c
834 
835  IF(k1 == k2) THEN
836  fpt = f_upper
837 
838  ELSE
839 
840  f0 = aw0(i,1)*field(n1,k2)+aw0(i,2)*field(n2,k2)+aw0(i,3)*field(n3,k2)
841  fx = awx(i,1)*field(n1,k2)+awx(i,2)*field(n2,k2)+awx(i,3)*field(n3,k2)
842  fy = awy(i,1)*field(n1,k2)+awy(i,2)*field(n2,k2)+awy(i,3)*field(n3,k2)
843  f_lower = f0 + fx*x0c + fy*y0c
844 
845  fpt = (alpha)*f_lower + (1.0_dp-alpha)*f_upper
846  END IF
847 
848 
849  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"END: INTERP_PNODAL_3D"
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target dzz1
Definition: mod_main.f90:1097
real(sp), dimension(:,:), allocatable, target awx
Definition: mod_main.f90:1333
real(sp), dimension(:,:), allocatable, target aw0
Definition: mod_main.f90:1335
real(sp), dimension(:,:), allocatable, target awy
Definition: mod_main.f90:1334
integer kb
Definition: mod_main.f90:64
integer kbm2
Definition: mod_main.f90:66
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
real(sp), dimension(:,:), allocatable, target zz1
Definition: mod_main.f90:1095
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
real(sp), dimension(:,:), allocatable, target z1
Definition: mod_main.f90:1094
integer kbm1
Definition: mod_main.f90:65

◆ interp_pnodal_3d_flt()

real(spa) function mod_utils::interp_pnodal::interp_pnodal_3d_flt ( real(spa), intent(in)  xloc,
real(spa), intent(in)  yloc,
real(spa), intent(in)  sigloc,
integer, intent(in)  lvls,
integer, intent(in)  i,
real(spa), dimension(:,:), intent(in), pointer  Field 
)

Definition at line 652 of file mod_utils.f90.

652  USE all_vars, only : aw0,awx,awy,nv,xc,yc,kbm2,kbm1,kb,z1,zz1,dz1,dzz1
653  IMPLICIT NONE
654  REAL(SPA) :: FPT
655  INTEGER, INTENT(IN) :: i, lvls ! Cell Number, Number of Levels in data (kb or kbm1)
656  REAL(SPA), INTENT(IN):: xloc,yloc,sigloc
657  REAL(SPA), POINTER, INTENT(IN) :: FIELD(:,:)
658 
659  REAL(SPA):: X0c,Y0c,F0,Fx,Fy,F_LOWER,F_UPPER, alpha,dsig
660  INTEGER :: n1,n2,n3,k1,k2,k
661  REAL(SPA) :: dx_sph
662 
663  IF(dbg_set(dbg_sbr)) WRITE(ipt,*) "START: INTERP_PNODAL_3D"
664 
665  !element location (i) and surrounding nodes (n1,n2,n3)
666  n1 = nv(i,1)
667  n2 = nv(i,2)
668  n3 = nv(i,3)
669 
670 !offset from element center
671 !---------------------------------------------------------------
672 
673  x0c = xloc - xc(i)
674  y0c = yloc - yc(i)
675 !---------------------------------------------------------------
676 !---------------------------------------------------------------
677 
678  ! Determine the layer in which the point resides
679  IF(lvls == kbm1) THEN
680 
681  !top
682  if(sigloc >= zz1(i,1))then
683  k1 = 1
684  k2 = 1
685  alpha = -1
686  elseif(sigloc > zz1(i,kbm1)) then!intermediate
687  do k=1,kbm2
688  if(sigloc < zz1(i,k) .and. sigloc >= zz1(i,k+1) )then
689  k1 = k
690  k2 = k+1
691  alpha = (zz1(i,k)-sigloc)/dzz1(i,k)
692  exit
693  endif
694  end do
695  else
696  ! TOP
697  k1 = kbm1
698  k2 = kbm1
699  alpha = -1
700  endif
701  ELSE IF(lvls == kb) THEN
702 
703  !top
704  if(sigloc >= z1(i,1))then
705  k1 = 1
706  k2 = 1
707  alpha = -1
708  elseif(sigloc > z1(i,kb)) then !intermediate
709  do k=1,kbm1
710  if(sigloc < z1(i,k) .and. sigloc >= z1(i,k+1) )then
711  k1 = k
712  k2 = k+1
713  alpha = (z1(i,k)-sigloc)/dz1(i,k)
714  endif
715  end do
716  else
717  !bottom
718  k1 = kbm1
719  k2 = kbm1
720  alpha = -1
721  endif
722 
723  ELSE
724  CALL fatal_error("INTERP_PNODAL_3D: Invalid number of levels passed",&
725  & "(Must be equal to either KB or KBM1")
726  END IF
727 
728  !linear interpolation of Field
729  f0 = aw0(i,1)*field(n1,k1)+aw0(i,2)*field(n2,k1)+aw0(i,3)*field(n3,k1)
730  fx = awx(i,1)*field(n1,k1)+awx(i,2)*field(n2,k1)+awx(i,3)*field(n3,k1)
731  fy = awy(i,1)*field(n1,k1)+awy(i,2)*field(n2,k1)+awy(i,3)*field(n3,k1)
732  f_upper = f0 + fx*x0c + fy*y0c
733 
734  IF(k1 == k2) THEN
735  fpt = f_upper
736 
737  ELSE
738 
739  f0 = aw0(i,1)*field(n1,k2)+aw0(i,2)*field(n2,k2)+aw0(i,3)*field(n3,k2)
740  fx = awx(i,1)*field(n1,k2)+awx(i,2)*field(n2,k2)+awx(i,3)*field(n3,k2)
741  fy = awy(i,1)*field(n1,k2)+awy(i,2)*field(n2,k2)+awy(i,3)*field(n3,k2)
742  f_lower = f0 + fx*x0c + fy*y0c
743 
744  fpt = (alpha)*f_lower + (1.0_spa-alpha)*f_upper
745  END IF
746 
747 
748  IF(dbg_set(dbg_sbr)) WRITE(ipt,*)"END: INTERP_PNODAL_3D"
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:,:), allocatable, target dzz1
Definition: mod_main.f90:1097
real(sp), dimension(:,:), allocatable, target awx
Definition: mod_main.f90:1333
real(sp), dimension(:,:), allocatable, target aw0
Definition: mod_main.f90:1335
real(sp), dimension(:,:), allocatable, target awy
Definition: mod_main.f90:1334
integer kb
Definition: mod_main.f90:64
integer kbm2
Definition: mod_main.f90:66
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
real(sp), dimension(:,:), allocatable, target zz1
Definition: mod_main.f90:1095
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:,:), allocatable, target dz1
Definition: mod_main.f90:1096
real(sp), dimension(:,:), allocatable, target z1
Definition: mod_main.f90:1094
integer kbm1
Definition: mod_main.f90:65

The documentation for this interface was generated from the following file: