ddx 0.6.8
Libary for domain-decomposition methods for polarizable continuum models
ddx_pcm.f90
Go to the documentation of this file.
1
11
13module ddx_pcm
14! Get ddx-operators
16implicit none
17
20
21contains
22
31subroutine pcm_energy(constants, state, esolv, ddx_error)
32 implicit none
33 type(ddx_constants_type), intent(in) :: constants
34 type(ddx_state_type), intent(in) :: state
35 type(ddx_error_type), intent(inout) :: ddx_error
36 real(dp), intent(out) :: esolv
37 real(dp), external :: ddot
38 ! dummy operation on unused interface arguments
39 if (ddx_error % flag .eq. 0) continue
40 esolv = pt5*ddot(constants % n, state % xs, 1, state % psi, 1)
41end subroutine pcm_energy
42
55subroutine pcm_setup(params, constants, workspace, state, phi_cav, psi, ddx_error)
56 implicit none
57 type(ddx_params_type), intent(in) :: params
58 type(ddx_constants_type), intent(in) :: constants
59 type(ddx_workspace_type), intent(inout) :: workspace
60 type(ddx_state_type), intent(inout) :: state
61 type(ddx_error_type), intent(inout) :: ddx_error
62 real(dp), intent(in) :: phi_cav(constants % ncav)
63 real(dp), intent(in) :: psi(constants % nbasis, params % nsph)
64
65 ! dummy operation on unused interface arguments
66 if (ddx_error % flag .eq. 0) continue
67
68 call cav_to_spherical(params, constants, workspace, phi_cav, &
69 & state % phi)
70 state % phi = - state % phi
71 state % phi_cav = phi_cav
72 state % psi = psi
73 state % rhs_done = .true.
74 state % adjoint_rhs_done = .true.
75end subroutine pcm_setup
76
85subroutine pcm_guess(params, constants, workspace, state, ddx_error)
86 implicit none
87 type(ddx_params_type), intent(in) :: params
88 type(ddx_constants_type), intent(in) :: constants
89 type(ddx_workspace_type), intent(inout) :: workspace
90 type(ddx_state_type), intent(inout) :: state
91 type(ddx_error_type), intent(inout) :: ddx_error
92
93 state % xs = zero
94 call prec_repsx(params, constants, workspace, state % phi, &
95 & state % phieps, ddx_error)
96
97end subroutine pcm_guess
98
107subroutine pcm_guess_adjoint(params, constants, workspace, state, ddx_error)
108 implicit none
109 type(ddx_params_type), intent(in) :: params
110 type(ddx_constants_type), intent(in) :: constants
111 type(ddx_workspace_type), intent(inout) :: workspace
112 type(ddx_state_type), intent(inout) :: state
113 type(ddx_error_type), intent(inout) :: ddx_error
114
115 state % y = zero
116 call ldm1x(params, constants, workspace, state % psi, state % s, ddx_error)
117
118end subroutine pcm_guess_adjoint
119
130subroutine pcm_solve(params, constants, workspace, state, tol, ddx_error)
131 implicit none
132 type(ddx_params_type), intent(in) :: params
133 type(ddx_constants_type), intent(in) :: constants
134 type(ddx_workspace_type), intent(inout) :: workspace
135 type(ddx_state_type), intent(inout) :: state
136 type(ddx_error_type), intent(inout) :: ddx_error
137 real(dp), intent(in) :: tol
138 ! local variables
139 real(dp) :: start_time, finish_time
140
141 call rinfx(params, constants, workspace, state % phi, state % phiinf, &
142 & ddx_error)
143
144 state % phieps_niter = params % maxiter
145 start_time = omp_get_wtime()
146 call jacobi_diis(params, constants, workspace, tol, state % phiinf, &
147 & state % phieps, state % phieps_niter, state % phieps_rel_diff, &
148 & repsx, prec_repsx, hnorm, ddx_error)
149 finish_time = omp_get_wtime()
150 state % phieps_time = finish_time - start_time
151
152 if (ddx_error % flag .ne. 0) then
153 call update_error(ddx_error, "ddpcm_solve: solver for ddPCM " // &
154 & "system did not converge, exiting")
155 return
156 end if
157
158 state % xs_niter = params % maxiter
159 start_time = omp_get_wtime()
160 call jacobi_diis(params, constants, workspace, tol, state % phieps, &
161 & state % xs, state % xs_niter, state % xs_rel_diff, lx, ldm1x, &
162 & hnorm, ddx_error)
163 finish_time = omp_get_wtime()
164 state % xs_time = finish_time - start_time
165
166 if (ddx_error % flag .ne. 0) then
167 call update_error(ddx_error, "ddpcm_solve: solver for ddCOSMO " // &
168 & "system did not converge, exiting")
169 return
170 end if
171
172end subroutine pcm_solve
173
184subroutine pcm_solve_adjoint(params, constants, workspace, state, tol, ddx_error)
185 implicit none
186 type(ddx_params_type), intent(in) :: params
187 type(ddx_constants_type), intent(in) :: constants
188 type(ddx_workspace_type), intent(inout) :: workspace
189 type(ddx_state_type), intent(inout) :: state
190 real(dp), intent(in) :: tol
191 type(ddx_error_type), intent(inout) :: ddx_error
192 ! local variables
193 real(dp) :: start_time, finish_time
194
195 state % s_niter = params % maxiter
196 start_time = omp_get_wtime()
197 call jacobi_diis(params, constants, workspace, tol, state % psi, &
198 & state % s, state % s_niter, state % s_rel_diff, lstarx, ldm1x, &
199 & hnorm, ddx_error)
200 finish_time = omp_get_wtime()
201 state % s_time = finish_time - start_time
202
203 if (ddx_error % flag .ne. 0) then
204 call update_error(ddx_error, "ddpcm_solve_adjoint: solver for ddCOSMO " // &
205 & "system did not converge, exiting")
206 return
207 end if
208
209 state % y_niter = params % maxiter
210 start_time = omp_get_wtime()
211 call jacobi_diis(params, constants, workspace, tol, state % s, state % y, &
212 & state % y_niter, state % y_rel_diff, repsstarx, prec_repsstarx, &
213 & hnorm, ddx_error)
214 finish_time = omp_get_wtime()
215 state % y_time = finish_time - start_time
216
217 if (ddx_error % flag .ne. 0) then
218 call update_error(ddx_error, "ddpcm_solve_adjoint: solver for ddPCM " // &
219 & "system did not converge, exiting")
220 return
221 end if
222
223 ! compute the real adjoint solution and store it in Q
224 state % q = state % s - fourpi/(params % eps - one)*state % y
225
226 call pcm_derivative_setup(params, constants, workspace, state)
227
228end subroutine pcm_solve_adjoint
229
241subroutine pcm_solvation_force_terms(params, constants, workspace, &
242 & state, e_cav, force, ddx_error)
243 implicit none
244 type(ddx_params_type), intent(in) :: params
245 type(ddx_constants_type), intent(in) :: constants
246 type(ddx_workspace_type), intent(inout) :: workspace
247 type(ddx_state_type), intent(inout) :: state
248 type(ddx_error_type), intent(inout) :: ddx_error
249 real(dp), intent(in) :: e_cav(3, constants % ncav)
250 real(dp), intent(out) :: force(3, params % nsph)
251
252 real(dp) :: start_time, finish_time
253 integer :: isph
254
255 ! dummy operation on unused interface arguments
256 if (ddx_error % flag .eq. 0) continue
257
258 start_time = omp_get_wtime()
259
260 call gradr(params, constants, workspace, state % g, state % ygrid, force)
261
262 do isph = 1, params % nsph
263 call contract_grad_l(params, constants, isph, state % xs, &
264 & state % sgrid, workspace % tmp_vylm(:, 1), &
265 & workspace % tmp_vdylm(:, :, 1), workspace % tmp_vplm(:, 1), &
266 & workspace % tmp_vcos(:, 1), workspace % tmp_vsin(:, 1), &
267 & force(:, isph))
268 call contract_grad_u(params, constants, isph, state % qgrid, &
269 & state % phi_grid, force(:, isph))
270 end do
271 force = - pt5 * force
272
273 call zeta_grad(params, constants, state, e_cav, force)
274
275 finish_time = omp_get_wtime()
276 state % force_time = finish_time - start_time
277
278end subroutine pcm_solvation_force_terms
279
288subroutine pcm_derivative_setup(params, constants, workspace, state)
289 type(ddx_params_type), intent(in) :: params
290 type(ddx_constants_type), intent(in) :: constants
291 type(ddx_workspace_type), intent(inout) :: workspace
292 type(ddx_state_type), intent(inout) :: state
293
294 real(dp), external :: ddot
295 integer :: icav, isph, igrid
296
297 ! dummy operation on unused interface arguments
298 if (allocated(workspace % tmp_pot)) continue
299
300 ! Get grid values of S and Y
301 call dgemm('T', 'N', params % ngrid, params % nsph, &
302 & constants % nbasis, one, constants % vgrid, constants % vgrid_nbasis, &
303 & state % s, constants % nbasis, zero, state % sgrid, params % ngrid)
304 call dgemm('T', 'N', params % ngrid, params % nsph, &
305 & constants % nbasis, one, constants % vgrid, constants % vgrid_nbasis, &
306 & state % y, constants % nbasis, zero, state % ygrid, params % ngrid)
307
308 ! Get the values of phi on the grid
309 call ddcav_to_grid_work(params % ngrid, params % nsph, constants % ncav, &
310 & constants % icav_ia, constants % icav_ja, state % phi_cav, &
311 & state % phi_grid)
312
313 state % g = - (state % phieps - state % phi)
314 state % qgrid = state % sgrid - fourpi/(params % eps - one)*state % ygrid
315
316 ! assemble the intermediate zeta: S weighted by U evaluated on the
317 ! exposed grid points.
318 icav = 0
319 do isph = 1, params % nsph
320 do igrid = 1, params % ngrid
321 if (constants % ui(igrid, isph) .ne. zero) then
322 icav = icav + 1
323 state % zeta(icav) = constants % wgrid(igrid) * &
324 & constants % ui(igrid, isph) * ddot(constants % nbasis, &
325 & constants % vgrid(1, igrid), 1, state % q(1, isph), 1)
326 end if
327 end do
328 end do
329
330end subroutine pcm_derivative_setup
331
332end module ddx_pcm
subroutine pcm_setup(params, constants, workspace, state, phi_cav, psi, ddx_error)
Given the potential at the cavity points, assemble the RHS for ddCOSMO or for ddPCM.
Definition: ddx_pcm.f90:56
subroutine pcm_guess(params, constants, workspace, state, ddx_error)
Do a guess for the primal ddPCM linear system.
Definition: ddx_pcm.f90:86
subroutine pcm_energy(constants, state, esolv, ddx_error)
Compute the ddPCM energy.
Definition: ddx_pcm.f90:32
subroutine pcm_solvation_force_terms(params, constants, workspace, state, e_cav, force, ddx_error)
Compute the solvation contribution to the ddPCM forces.
Definition: ddx_pcm.f90:243
subroutine pcm_solve(params, constants, workspace, state, tol, ddx_error)
Solve the ddPCM primal linear system.
Definition: ddx_pcm.f90:131
subroutine pcm_solve_adjoint(params, constants, workspace, state, tol, ddx_error)
Solve the ddPCM adjpint linear system.
Definition: ddx_pcm.f90:185
subroutine pcm_guess_adjoint(params, constants, workspace, state, ddx_error)
Do a guess for the adjoint ddPCM linear system.
Definition: ddx_pcm.f90:108
Operators shared among ddX methods.
subroutine repsx(params, constants, workspace, x, y, ddx_error)
Apply operator to spherical harmonics.
subroutine prec_repsstarx(params, constants, workspace, x, y, ddx_error)
Apply preconditioner for ddPCM adjoint linear system.
subroutine prec_repsx(params, constants, workspace, x, y, ddx_error)
Apply preconditioner for ddPCM primal linear system.
subroutine lx(params, constants, workspace, x, y, ddx_error)
Single layer operator matvec.
subroutine rinfx(params, constants, workspace, x, y, ddx_error)
Apply operator to spherical harmonics.
subroutine ldm1x(params, constants, workspace, x, y, ddx_error)
Diagonal preconditioning for Lx operator.
subroutine repsstarx(params, constants, workspace, x, y, ddx_error)
Apply operator to spherical harmonics.
subroutine lstarx(params, constants, workspace, x, y, ddx_error)
Adjoint single layer operator matvec.
High-level subroutines for ddpcm.
Definition: ddx_pcm.f90:13
subroutine pcm_derivative_setup(params, constants, workspace, state)
This routines precomputes the intermediates to be used in the evaluation of ddCOSMO analytical deriva...
Definition: ddx_pcm.f90:289