ddx 0.6.0
Libary for domain-decomposition methods for polarizable continuum models
ddx_parameters.f90
Go to the documentation of this file.
1
13
16! Include compile-time definitions
18use ddx_errors
19! Enable OpenMP
20use omp_lib
21implicit none
22
23
27 integer :: model
29 integer :: force
31 real(dp) :: eps
33 real(dp) :: kappa
35 real(dp) :: eta
38 real(dp) :: se
40 integer :: lmax
42 integer :: ngrid
44 integer :: maxiter
47 integer :: jacobi_ndiis
49 integer :: fmm
52 integer :: pm
55 integer :: pl
58 integer :: nproc
60 integer :: nsph
62 real(dp), allocatable :: csph(:, :)
64 real(dp), allocatable :: rsph(:)
66 real(dp) :: epsp = 1.0_dp
68 integer :: matvecmem
71 character(len=255) :: output_filename
73 integer :: len_output_filename
75 logical :: verbose
77 integer :: iunit
78end type ddx_params_type
79
80contains
81
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)
123 !! Inputs
124 ! Model to use 1 for COSMO, 2 for PCM, 3 for LPB.
125 integer, intent(in) :: model
126 ! Whether computing analytical forces will be required (1) or not (0).
127 integer, intent(in) :: force
128 ! Relative dielectric permittivity.
129 real(dp), intent(in) :: eps
130 ! Debye-H\"{u}ckel parameter. Referenced only in LPB model (model=3)
131 real(dp), intent(in) :: kappa
132 ! Regularization parameter.
133 real(dp), intent(in) :: eta
134 ! Shift of the regularization. -1 for interior, 0 for centered and
135 ! 1 for outer regularization.
136 real(dp), intent(in) :: se
137 ! Maximal degree of modeling spherical harmonics.
138 integer, intent(in) :: lmax
139 ! Number of Lebedev grid points on each sphere.
140 integer, intent(in) :: ngrid
141 ! handling of sparse matrix. 1 for precomputing them and keeping them in
142 ! memory, 0 for assembling the mvps on-the-fly.
143 integer, intent(in) :: matvecmem
144 ! Maximum number of iterations for the iterative solver.
145 integer, intent(in) :: maxiter
146 ! Number of extrapolation points for Jacobi/DIIS solver. Referenced only
147 ! if Jacobi solver is used.
148 integer, intent(in) :: jacobi_ndiis
149 ! Enable (1) or disable (0) use of FMM techniques.
150 integer, intent(in) :: fmm
151 ! Maximal degree of spherical harmonics for a multipole expansion.
152 !
153 ! If this value is -1 then no far-field FMM interactions are performed.
154 integer, intent(in) :: pm
155 ! Maximal degree of spherical harmonics for a local expansion.
156 !
157 ! If this value is -1 then no far-field FMM interactions are performed.
158 integer, intent(in) :: pl
159 ! Number of OpenMP threads to be used.
160 integer, intent(in) :: nproc
161 ! Number of atoms in the molecule.
162 integer, intent(in) :: nsph
163 ! Centers of atoms of a dimension (3, nsph).
164 real(dp), intent(in) :: csph(3, nsph)
165 ! Array of radii of atoms of a dimension (nsph).
166 real(dp), intent(in) :: rsph(nsph)
167 ! log file name
168 character(len=255) :: output_filename
169 !! Outputs
170 type(ddx_params_type), intent(out) :: params
171 type(ddx_error_type), intent(inout) :: ddx_error
172 !! Local variables
173 integer :: igrid, i, info
174 !! The code
175 if (ddx_error % flag .ne. 0) then
176 call update_error(ddx_error, "params_init received input in ddx_error " // &
177 & " state, exiting")
178 end if
179 ! parse the log file name
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))
184 else
185 params % output_filename = ''
186 params % verbose = .false.
187 params % len_output_filename = 0
188 end if
189 ! Model, 1=COSMO, 2=PCM, 3=LPB
190 if ((model .lt. 1) .or. (model .gt. 3)) then
191 call update_error(ddx_error, "params_init: invalid value of `model`")
192 end if
193 params % model = model
194 ! Check if forces are needed
195 if ((force .lt. 0) .or. (force .gt. 1)) then
196 call update_error(ddx_error, "params_init: invalid value of `force`")
197 end if
198 params % force = force
199 ! Relative dielectric permittivity
200 if (eps .le. one) then
201 call update_error(ddx_error, "params_init: invalid value of `eps`")
202 end if
203 params % eps = eps
204 ! Debye-H\"{u}ckel parameter (only used in ddLPB)
205 if ((model .eq. 3) .and. (kappa .le. zero)) then
206 call update_error(ddx_error, "params_init: invalid value of `kappa`")
207 end if
208 params % kappa = kappa
209 ! Regularization parameter
210 if ((eta .lt. zero) .or. (eta .gt. one)) then
211 call update_error(ddx_error, "params_init: invalid value of `eta`")
212 end if
213 params % eta = eta
214 ! Shift of a regularization
215 if ((se .lt. -one) .or. (se .gt. one)) then
216 call update_error(ddx_error, "params_init: invalid value of `se`")
217 end if
218 params % se = se
219 ! Degree of modeling spherical harmonics
220 if (lmax .lt. 0) then
221 call update_error(ddx_error, "params_init: invalid value of `lmax`")
222 end if
223 params % lmax = lmax
224 ! Check number of Lebedev grid points
225 igrid = 0
226 do i = 1, nllg
227 if (ng0(i) .eq. ngrid) then
228 igrid = i
229 exit
230 end if
231 end do
232 if (igrid .eq. 0) then
233 call update_error(ddx_error, "params_init: Unsupported value of `ngrid`")
234 end if
235 params % ngrid = ngrid
236 ! Maximum number of iterations
237 if (maxiter .le. 0) then
238 call update_error(ddx_error, "params_init: invalid value of `maxiter`")
239 end if
240 params % maxiter = maxiter
241 ! Number of Jacobi DIIS extrapolation points (ndiis=25 works)
242 if (jacobi_ndiis .lt. 0) then
243 call update_error(ddx_error, "params_init: invalid value of `jacobi_ndiis`")
244 end if
245 params % jacobi_ndiis = jacobi_ndiis
246 ! Check if FMM-acceleration is needed
247 if ((fmm .lt. 0) .or. (fmm .gt. 1)) then
248 call update_error(ddx_error, "params_init: invalid value of `fmm`")
249 end if
250 params % fmm = fmm
251 ! Set FMM parameters if FMM is needed
252 if (fmm .eq. 1) then
253 ! Maximal degree of multipole spherical harmonics. Value -1 means no
254 ! far-field interactions are to be computed, only near-field
255 ! interactions are taken into account.
256 if (pm .lt. -1) then
257 call update_error(ddx_error, "params_init: invalid value of `pm`")
258 end if
259 ! Maximal degree of local spherical harmonics. Value -1 means no
260 ! far-field interactions are to be computed, only near-field
261 ! interactions are taken into account.
262 if (pl .lt. -1) then
263 call update_error(ddx_error, "params_init: invalid value of `pl`")
264 end if
265 ! If far-field interactions are to be ignored
266 if ((pl .eq. -1) .or. (pm .eq. -1)) then
267 params % pm = -1
268 params % pl = -1
269 ! If far-field interactions are to be taken into account
270 else
271 params % pm = pm
272 params % pl = pl
273 end if
274 else
275 ! These values are ignored if fmm flag is 0
276 params % pm = -2
277 params % pl = -2
278 end if
279 ! Number of OpenMP threads to be used
280 ! available.
281 if (nproc .lt. 0) then
282 call update_error(ddx_error, "params_init: invalid value of `nproc`")
283 params % nproc = 1
284 else if (nproc .eq. 0) then
285 params % nproc = 1
286 else
287 params % nproc = nproc
288 end if
289 call omp_set_num_threads(params % nproc)
290 ! Number of atoms
291 if (nsph .le. 0) then
292 call update_error(ddx_error, "params_init: invalid value of `nsph`")
293 end if
294 params % nsph = 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")
299 return
300 end if
301 params % csph = csph
302 params % rsph = rsph
303 if (matvecmem.eq.0 .or. matvecmem.eq.1) then
304 params % matvecmem = matvecmem
305 else
306 call update_error(ddx_error, "params_init: invalid value of `matvecmem`")
307 end if
308
309 if (ddx_error % flag .ne. 0) return
310
311 ! init log
312 call init_printing(params, ddx_error)
313 if (ddx_error % flag .ne. 0) then
314 call update_error(ddx_error, "init_printing returned an error, exiting")
315 return
316 end if
317end subroutine params_init
318
324 integer, intent(inout) :: ngrid
325 integer :: igrid, i, inear, jnear
326 ! Get nearest number of Lebedev grid points
327 igrid = 0
328 inear = 100000
329 do i = 1, nllg
330 jnear = abs(ng0(i) - ngrid)
331 if (jnear .lt. inear) then
332 inear = jnear
333 igrid = i
334 end if
335 end do
336 ngrid = ng0(igrid)
337end subroutine
338
344subroutine params_free(params, ddx_error)
345 implicit none
346 type(ddx_params_type), intent(inout) :: params
347 type(ddx_error_type), intent(inout) :: ddx_error
348 integer :: istat
349
350 call finalize_printing(params)
351
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")
356 end if
357 end if
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")
362 end if
363 end if
364end subroutine params_free
365
367subroutine init_printing(params, ddx_error)
368 implicit none
369 type(ddx_params_type), intent(inout) :: params
370 type(ddx_error_type), intent(inout) :: ddx_error
371 logical :: exists
372 if (.not.params % verbose) return
373 inquire(file=params % output_filename(1:params % len_output_filename), &
374 & exist=exists)
375 if (exists) then
376 call update_error(ddx_error, 'Log file already present')
377 return
378 else
379 params % iunit = 100
380 open(params % iunit, &
381 & file=params % output_filename(1: params % len_output_filename), &
382 & form='formatted')
383 end if
384end subroutine init_printing
385
387subroutine finalize_printing(params)
388 implicit none
389 type(ddx_params_type), intent(inout) :: params
390 if (.not.params % verbose) return
391 close(params % iunit)
392 params % verbose = .false.
393 params % output_filename = ''
394 params % len_output_filename = 0
395 params % iunit = 0
396end subroutine finalize_printing
397
398end module ddx_parameters
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
Definition: ddx_errors.f90:6
Type to check and store user input parameters.