ObserverConstruct Subroutine

public subroutine ObserverConstruct(self, sourceZone, mesh, ID, solution_file, FirstCall, interpolate, totalNumberOfFaces, elementSide)

Uses

    • ParamfileRegions
    • FileReadingUtilities

Type Bound

ObserverClass

Arguments

Type IntentOptional Attributes Name
class(ObserverClass) :: self
class(Zone_t), intent(in) :: sourceZone
class(HexMesh), intent(in) :: mesh
integer, intent(in) :: ID
character(len=*), intent(in) :: solution_file
logical, intent(in) :: FirstCall
logical, intent(in) :: interpolate
integer, intent(in) :: totalNumberOfFaces
integer, intent(in), dimension(:) :: elementSide

Source Code

   Subroutine ObserverConstruct(self, sourceZone, mesh, ID, solution_file, FirstCall, interpolate, totalNumberOfFaces, elementSide)

!        *****************************************************************************
!              This subroutine initializes the observer similar to a monitor. The following
!           data is obtained from the case file:
!              -> Name: The observer name (10 characters maximum)
!              -> x: The observer position
!        *****************************************************************************

       use ParamfileRegions
       use FileReadingUtilities, only: getRealArrayFromString
       implicit none

       class(ObserverClass)                                 :: self
       class(Zone_t), intent(in)                            :: sourceZone
       class(HexMesh), intent(in)                           :: mesh
       integer, intent(in)                                  :: ID, totalNumberOfFaces
       character(len=*), intent(in)                         :: solution_file
       logical, intent(in)                                  :: FirstCall, interpolate
       integer, dimension(:), intent(in)                    :: elementSide

       ! local variables
       character(len=STR_LEN_OBSERVER)  :: in_label
       character(len=STR_LEN_OBSERVER)  :: fileName
       character(len=STR_LEN_OBSERVER)  :: paramFile
       character(len=STR_LEN_OBSERVER)  :: coordinates
       integer                          :: fID
       integer                          :: MeshFaceID, zoneFaceID
       ! integer                          :: elementSide
!
!      Get observer ID
!      --------------
       self % ID = ID
!
!      Search for the parameters in the case file
!      ------------------------------------------
       write(in_label , '(A,I0)') "#define acoustic observer " , self % ID

       call get_command_argument(1, paramFile)
       call readValueInRegion(trim ( paramFile), "name",   self % observerName, in_label, "# end" ) 
       call readValueInRegion(trim(paramFile), "position", coordinates        , in_label, "# end" )

!      Get the coordinates
!      -------------------
       self % x = getRealArrayFromString(coordinates)

!     Enable the observer
!     ------------------
      self % active = .true.
      allocate ( self % Pac(OB_BUFFER_SIZE,3) )

!     Get source information
!     ------------------
      self % numberOfFaces = sourceZone % no_of_faces

!     Construct each pair observer-source
!     ------------------
      allocate( self % sourcePair(self % numberOfFaces) )
!     Loop the zone to get faces
      do zoneFaceID = 1, self % numberOfFaces
!         Face global ID
          MeshFaceID = sourceZone % faces(zoneFaceID)
          call self % sourcePair(zoneFaceID) % construct(self % x, mesh % faces(MeshFaceID), MeshFaceID, FirstCall, elementSide(zoneFaceID))
      end do  

!     Allocate variables for interpolation
!     -------------------------------------------------
      if (interpolate) then
          do zoneFaceID = 1, self % numberOfFaces
              call self % sourcePair(zoneFaceID) % allocPacc(OB_BUFFER_SIZE)
          end do
      end if 

!     Set the average time delay of the observer
!     -------------------------------------------------
      call self % updateTdelay(totalNumberOfFaces)

!     Prepare the file in which the observer is exported
!     -------------------------------------------------
      write( self % fileName , '(A,A,A,A)') trim(solution_file) , "." , trim(self % observerName) , ".observer"
!
!     Create file
!     -----------
      if (FirstCall) then
         open ( newunit = fID , file = trim(self % fileName) , status = "unknown" , action = "write" ) 

!        Write the file headers
!        ----------------------
         write( fID , '(A20,A  )') "Observer name:      ", trim(self % observerName)
         write( fID , '(A20,ES24.10)') "x coordinate: ", self % x(1)
         write( fID , '(A20,ES24.10)') "y coordinate: ", self % x(2)
         write( fID , '(A20,ES24.10)') "z coordinate: ", self % x(3)

         write( fID , * )
         write( fID , '(A10,5(2X,A24))' ) "Iteration" , "Time" , "Observer_Time", "P'T", "P'L", "P'"

         close ( fID )
      end if

   End Subroutine ObserverConstruct