47 integer :: jacobi_ndiis
62 real(
dp),
allocatable :: csph(:, :)
64 real(
dp),
allocatable :: rsph(:)
66 real(
dp) :: epsp = 1.0_dp
71 character(len=255) :: output_filename
73 integer :: len_output_filename
120subroutine params_init(model, force, eps, kappa, eta, se, lmax, ngrid, &
121 & matvecmem, maxiter, jacobi_ndiis, fmm, pm, pl, nproc, nsph, &
122 & csph, rsph, output_filename, params, ddx_error)
125 integer,
intent(in) :: model
127 integer,
intent(in) :: force
129 real(dp),
intent(in) :: eps
131 real(dp),
intent(in) :: kappa
133 real(dp),
intent(in) :: eta
136 real(dp),
intent(in) :: se
138 integer,
intent(in) :: lmax
140 integer,
intent(in) :: ngrid
143 integer,
intent(in) :: matvecmem
145 integer,
intent(in) :: maxiter
148 integer,
intent(in) :: jacobi_ndiis
150 integer,
intent(in) :: fmm
154 integer,
intent(in) :: pm
158 integer,
intent(in) :: pl
160 integer,
intent(in) :: nproc
162 integer,
intent(in) :: nsph
164 real(dp),
intent(in) :: csph(3, nsph)
166 real(dp),
intent(in) :: rsph(nsph)
168 character(len=255) :: output_filename
173 integer :: igrid, i, info
175 if (ddx_error % flag .ne. 0)
then
176 call update_error(ddx_error,
"params_init received input in ddx_error " // &
180 if (len(trim(output_filename)) .ne. 0)
then
181 params % output_filename = output_filename
182 params % verbose = .true.
183 params % len_output_filename = len(trim(output_filename))
185 params % output_filename =
''
186 params % verbose = .false.
187 params % len_output_filename = 0
190 if ((model .lt. 1) .or. (model .gt. 3))
then
191 call update_error(ddx_error,
"params_init: invalid value of `model`")
193 params % model = model
195 if ((force .lt. 0) .or. (force .gt. 1))
then
196 call update_error(ddx_error,
"params_init: invalid value of `force`")
198 params % force = force
200 if (eps .le. one)
then
201 call update_error(ddx_error,
"params_init: invalid value of `eps`")
205 if ((model .eq. 3) .and. (kappa .le. zero))
then
206 call update_error(ddx_error,
"params_init: invalid value of `kappa`")
208 params % kappa = kappa
210 if ((eta .lt. zero) .or. (eta .gt. one))
then
211 call update_error(ddx_error,
"params_init: invalid value of `eta`")
215 if ((se .lt. -one) .or. (se .gt. one))
then
216 call update_error(ddx_error,
"params_init: invalid value of `se`")
220 if (lmax .lt. 0)
then
221 call update_error(ddx_error,
"params_init: invalid value of `lmax`")
227 if (
ng0(i) .eq. ngrid)
then
232 if (igrid .eq. 0)
then
233 call update_error(ddx_error,
"params_init: Unsupported value of `ngrid`")
235 params % ngrid = ngrid
237 if (maxiter .le. 0)
then
238 call update_error(ddx_error,
"params_init: invalid value of `maxiter`")
240 params % maxiter = maxiter
242 if (jacobi_ndiis .lt. 0)
then
243 call update_error(ddx_error,
"params_init: invalid value of `jacobi_ndiis`")
245 params % jacobi_ndiis = jacobi_ndiis
247 if ((fmm .lt. 0) .or. (fmm .gt. 1))
then
248 call update_error(ddx_error,
"params_init: invalid value of `fmm`")
257 call update_error(ddx_error,
"params_init: invalid value of `pm`")
263 call update_error(ddx_error,
"params_init: invalid value of `pl`")
266 if ((pl .eq. -1) .or. (pm .eq. -1))
then
281 if (nproc .lt. 0)
then
282 call update_error(ddx_error,
"params_init: invalid value of `nproc`")
284 else if (nproc .eq. 0)
then
287 params % nproc = nproc
289 call omp_set_num_threads(params % nproc)
291 if (nsph .le. 0)
then
292 call update_error(ddx_error,
"params_init: invalid value of `nsph`")
295 allocate(params % csph(3, nsph), params % rsph(nsph), stat=info)
296 if (info .ne. 0)
then
297 call update_error(ddx_error,
"params_init: `csph` and `rsph` " // &
298 &
"allocations failed")
303 if (matvecmem.eq.0 .or. matvecmem.eq.1)
then
304 params % matvecmem = matvecmem
306 call update_error(ddx_error,
"params_init: invalid value of `matvecmem`")
309 if (ddx_error % flag .ne. 0)
return
313 if (ddx_error % flag .ne. 0)
then
314 call update_error(ddx_error,
"init_printing returned an error, exiting")
324 integer,
intent(inout) :: ngrid
325 integer :: igrid, i, inear, jnear
330 jnear = abs(
ng0(i) - ngrid)
331 if (jnear .lt. inear)
then
352 if (
allocated(params % csph))
then
353 deallocate(params % csph, stat=istat)
354 if (istat .ne. 0)
then
355 call update_error(ddx_error,
"params_free: `csph` deallocation failed")
358 if (
allocated(params % rsph))
then
359 deallocate(params % rsph, stat=istat)
360 if (istat .ne. 0)
then
361 call update_error(ddx_error,
"params_free: `rsph` deallocation failed")
372 if (.not.params % verbose)
return
373 inquire(file=params % output_filename(1:params % len_output_filename), &
376 call update_error(ddx_error,
'Log file already present')
380 open(params % iunit, &
381 & file=params % output_filename(1: params % len_output_filename), &
390 if (.not.params % verbose)
return
391 close(params % iunit)
392 params % verbose = .false.
393 params % output_filename =
''
394 params % len_output_filename = 0
Compile-time constants and definitions.
integer, dimension(nllg), parameter ng0
Number of grid points of each Lebedev grid.
integer, parameter dp
Kind of double precision.
integer, parameter nllg
Number of supported Lebedev grids.
Module to treat properly user input parameters.
subroutine closest_supported_lebedev_grid(ngrid)
Adjust a guess for the number of Lebedev grid points.
subroutine init_printing(params, ddx_error)
Open the log file.
subroutine finalize_printing(params)
Close the log file.
subroutine params_free(params, ddx_error)
Deallocate the parameter object.
subroutine params_init(model, force, eps, kappa, eta, se, lmax, ngrid, matvecmem, maxiter, jacobi_ndiis, fmm, pm, pl, nproc, nsph, csph, rsph, output_filename, params, ddx_error)
Initialize and check input parameters.
ddX type for containing ddx_error information
Type to check and store user input parameters.