static_data_field_solver_s.f90 Source File


Contents


Source Code

submodule (static_data_m) static_data_field_solver_s
    !! Handles static data related with field_solver
    use params_hsolver_m, only : read_params_hsolver, write_params_hsolver, &
        params_hsolver, select_hsolver
    use device_handling_m, only: impose_default_device_affinity
    use descriptors_m, only : BND_TYPE_DIRICHLET_ZERO
    use helmholtz_solver_factory_m, only : helmholtz_solver_factory
    implicit none

contains

    module subroutine static_field_solver_init(par_filepath)
        character(len=*), intent(in), optional :: par_filepath

        character(len=:), allocatable :: fpath
        real(GP), allocatable, dimension(:) :: co, lambda, xi
        ! Dummy coefficients

        if (is_rank_info_writer) then
            write(get_stdout(),*)''
            write(get_stdout(),*)repeat('-',80)
            write(get_stdout(),*)'Initializing static data: field_solver'
        endif

        if (.not.equi_is_initialized) then
            call handle_error('equi must be initialized at field_solver initialization', &
                              GRILLIX_ERR_STATIC_DATA, __LINE__, __FILE__)
        endif
        if (.not.multigrid_is_initialized) then
            call handle_error('multigrid must be initialized at field_solver initialization', &
                              GRILLIX_ERR_STATIC_DATA, __LINE__, __FILE__)
        endif

        if (present(par_filepath)) then
            fpath = par_filepath
        else
            fpath = static_data_parpath
        endif

        ! Read parameters for field solver (PARALLAX functionality)
        call read_params_hsolver(fpath)
        if (is_rank_info_writer) call write_params_hsolver()

        ! Choose backend device
        call impose_default_device_affinity(select_hsolver)

        ! Create field solver for canonical mesh
        allocate(co(mesh_cano%get_n_points()), source=1.0_GP)
        allocate(lambda(mesh_cano%get_n_points_inner()), source=1.0_GP)
        allocate(xi(mesh_cano%get_n_points_inner()), source=1.0_GP)
        call helmholtz_solver_factory(multigrid_cano, &
                                      BND_TYPE_DIRICHLET_ZERO, &
                                      BND_TYPE_DIRICHLET_ZERO, &
                                      BND_TYPE_DIRICHLET_ZERO, &
                                      BND_TYPE_DIRICHLET_ZERO, &
                                      co, lambda, xi, &
                                      params_hsolver, select_hsolver, field_solver_cano)
        deallocate(xi)
        deallocate(lambda)
        deallocate(co)

        ! Create field solver for staggered mesh
        if (equi%is_axisymmetric()) then
            field_solver_stag => field_solver_cano
        else
            allocate(co(mesh_stag%get_n_points()),source=1.0_GP)
            allocate(lambda(mesh_stag%get_n_points_inner()),source=1.0_GP)
            allocate(xi(mesh_stag%get_n_points_inner()),source=1.0_GP)
            call helmholtz_solver_factory(multigrid_stag, &
                                          BND_TYPE_DIRICHLET_ZERO, &
                                          BND_TYPE_DIRICHLET_ZERO, &
                                          BND_TYPE_DIRICHLET_ZERO, &
                                          BND_TYPE_DIRICHLET_ZERO, &
                                          co, lambda, xi, &
                                          params_hsolver, select_hsolver, field_solver_stag_mem)
            deallocate(xi)
            deallocate(lambda)
            deallocate(co)

            field_solver_stag => field_solver_stag_mem
        endif

        field_solver_is_initialized = .true.

        if (is_rank_info_writer) then
            write(get_stdout(),*)repeat('-',80)
            write(get_stdout(),*)''
        endif

    end subroutine

end submodule