FWHConstruct Subroutine

public subroutine FWHConstruct(self, mesh, controlVariables)

Uses

    • FileReadingUtilities
    • mpi
    • FTValueDictionaryClass
    • mainKeywordsModule
    • Headers
    • Utilities
    • FWHDefinitions

Type Bound

FWHClass

Arguments

Type IntentOptional Attributes Name
class(FWHClass) :: self
class(HexMesh), intent(in) :: mesh
class(FTValueDictionary), intent(in) :: controlVariables

Source Code

    Subroutine FWHConstruct(self, mesh, controlVariables)

        use FTValueDictionaryClass
        use mainKeywordsModule
        use FileReadingUtilities, only: getCharArrayFromString
        use FWHDefinitions,       only: getMeanStreamValues
        use Headers
        use Utilities,            only: toLower
#ifdef _HAS_MPI_
        use mpi
#endif
        implicit none

        class(FWHClass)                                     :: self
        class(HexMesh), intent(in)                          :: mesh
        class(FTValueDictionary), intent(in)                :: controlVariables

!       ---------------
!       Local variables
!       ---------------
!
        integer                                             :: fID , io
        integer                                             :: i
        character(len=STR_LEN_OBSERVER)                     :: line
        character(len=STR_LEN_OBSERVER)                     :: solution_file
        integer                                             :: no_of_zones, no_of_face_i, ierr, no_of_faces
        logical, save                                       :: FirstCall = .TRUE.

!        look if the acoustic analogy calculations are set to be computed
!        --------------------------------
        !TODO read acoustic analogy type and return if is not defined, check for FWH if is defined and not FWH stop and send error
        if (.not. controlVariables % containsKey("acoustic analogy")) then
            self % isActive = .FALSE.
            ! print *, "FWH not activated"
            return
        end if

!       check the that sourceZone is FWH
!       ----------------------------------
        if (surfacesMesh % surfaceTypes(FWH_POSITION) .ne. SURFACE_TYPE_FWH) then
            self % isActive = .FALSE.
            print *, "FWH surface not found, the FWH routines will not deactivated"
            return
        end if 

!       Setup the buffer
!       ----------------
        if (controlVariables % containsKey("observers flush interval") ) then
           OB_BUFFER_SIZE = controlVariables % integerValueForKey("observers flush interval")
        end if

        self % isActive = .TRUE.
        allocate( self % t(OB_BUFFER_SIZE), self % iter(OB_BUFFER_SIZE) )

!       Get the general configuration of control file
!       First get the surface as a zone
!       -------------------------------
        self % isSolid   = .not. controlVariables % logicalValueForKey("acoustic analogy permeable")

!       Get the solution file name
!       --------------------------
        solution_file = controlVariables % stringValueForKey( solutionFileNameKey, requestedLength = STR_LEN_OBSERVER )
!
!       Remove the *.hsol termination
!       -----------------------------
        solution_file = trim(getFileName(solution_file))
        self % solution_file = trim(solution_file)

!       Search in case file for observers
!       ---------------------------------------------------------------------
        if (mesh % child) then ! Return doing nothing if this is a child mesh
           self % numberOfObservers = 0
        else
           self % numberOfObservers = getNoOfObservers()
        end if

!       Set interpolate attribute as TRUE by default
        ! todo: read from constrol variables
        self % interpolate = .TRUE.
        ! self % interpolate = .FALSE.

!       Initialize observers
!       --------------------
        call getMeanStreamValues()
        no_of_faces = surfacesMesh % totalFaces(SURFACE_TYPE_FWH)
        allocate( self % observers(self % numberOfObservers) )
        do i = 1, self % numberOfObservers
            call self % observers(i) % construct(surfacesMesh % zones(SURFACE_TYPE_FWH) , mesh, i, self % solution_file, FirstCall, &
                                                 self % interpolate, no_of_faces, surfacesMesh % elementSide(:,1))
        end do 

        self % bufferLine = 0
        self % firstWrite = .FALSE.
        
        FirstCall = .FALSE.

!        Describe the zones
!        ------------------
         if ( .not. MPI_Process % isRoot ) return
         call Subsection_Header("Fictitious FWH zone")
         write(STD_OUT,'(30X,A,A28,I0)') "->", "Number of faces: ", no_of_faces
         write(STD_OUT,'(30X,A,A28,I0)') "->", "Number of observers: ", self % numberOfObservers
         write(STD_OUT,'(30X,A,A28,I0)') "->", "Number of integrals: ", self % numberOfObservers * no_of_faces
         write(STD_OUT,'(30X,A,A28,L1)') "->", "Save zone solution: ", controlVariables % containsKey("acoustic save timestep")

    End Subroutine FWHConstruct