34 type(ddx_constants_type),
intent(in) :: constants
35 type(ddx_state_type),
intent(in) :: state
36 type(ddx_error_type),
intent(inout) :: ddx_error
37 real(dp),
intent(out) :: esolv
38 real(dp),
external :: ddot
41 if (ddx_error % flag .eq. 0)
continue
43 esolv = pt5*ddot(constants % n, state % xs, 1, state % psi, 1)
58subroutine cosmo_setup(params, constants, workspace, state, phi_cav, psi, ddx_error)
60 type(ddx_params_type),
intent(in) :: params
61 type(ddx_constants_type),
intent(in) :: constants
62 type(ddx_workspace_type),
intent(inout) :: workspace
63 type(ddx_state_type),
intent(inout) :: state
64 type(ddx_error_type),
intent(inout) :: ddx_error
65 real(dp),
intent(in) :: phi_cav(constants % ncav)
66 real(dp),
intent(in) :: psi(constants % nbasis, params % nsph)
69 if (ddx_error % flag .eq. 0)
continue
71 call cav_to_spherical(params, constants, workspace, phi_cav, &
73 state % phi = - state % phi
74 state % phi_cav = phi_cav
76 state % rhs_done = .true.
77 state % adjoint_rhs_done = .true.
89subroutine cosmo_guess(params, constants, workspace, state, ddx_error)
91 type(ddx_params_type),
intent(in) :: params
92 type(ddx_constants_type),
intent(in) :: constants
93 type(ddx_workspace_type),
intent(inout) :: workspace
94 type(ddx_state_type),
intent(inout) :: state
95 type(ddx_error_type),
intent(inout) :: ddx_error
98 call ldm1x(params, constants, workspace, state % phi, state % xs, ddx_error)
113 type(ddx_params_type),
intent(in) :: params
114 type(ddx_constants_type),
intent(in) :: constants
115 type(ddx_workspace_type),
intent(inout) :: workspace
116 type(ddx_state_type),
intent(inout) :: state
117 type(ddx_error_type),
intent(inout) :: ddx_error
120 call ldm1x(params, constants, workspace, state % psi, state % s, ddx_error)
134subroutine cosmo_solve(params, constants, workspace, state, tol, ddx_error)
136 type(ddx_params_type),
intent(in) :: params
137 type(ddx_constants_type),
intent(in) :: constants
138 type(ddx_workspace_type),
intent(inout) :: workspace
139 type(ddx_state_type),
intent(inout) :: state
140 real(dp),
intent(in) :: tol
141 type(ddx_error_type),
intent(inout) :: ddx_error
143 real(dp) :: start_time, finish_time
145 state % xs_niter = params % maxiter
146 start_time = omp_get_wtime()
147 call jacobi_diis(params, constants, workspace, tol, state % phi, &
148 & state % xs, state % xs_niter, state % xs_rel_diff,
lx,
ldm1x, &
150 finish_time = omp_get_wtime()
151 state % xs_time = finish_time - start_time
168 type(ddx_params_type),
intent(in) :: params
169 type(ddx_constants_type),
intent(in) :: constants
170 type(ddx_workspace_type),
intent(inout) :: workspace
171 type(ddx_state_type),
intent(inout) :: state
172 real(dp),
intent(in) :: tol
173 type(ddx_error_type),
intent(inout) :: ddx_error
175 real(dp) :: start_time, finish_time
177 state % s_niter = params % maxiter
178 start_time = omp_get_wtime()
179 call jacobi_diis(params, constants, workspace, tol, state % psi, &
180 & state % s, state % s_niter, state % s_rel_diff,
lstarx,
ldm1x, &
182 finish_time = omp_get_wtime()
183 state % s_time = finish_time - start_time
185 state % q = state % s
204 & state, e_cav, force, ddx_error)
206 type(ddx_params_type),
intent(in) :: params
207 type(ddx_constants_type),
intent(in) :: constants
208 type(ddx_workspace_type),
intent(inout) :: workspace
209 type(ddx_state_type),
intent(inout) :: state
210 type(ddx_error_type),
intent(inout) :: ddx_error
211 real(dp),
intent(in) :: e_cav(3, constants % ncav)
212 real(dp),
intent(inout) :: force(3, params % nsph)
214 real(dp) :: start_time, finish_time
218 if (ddx_error % flag .eq. 0)
continue
220 start_time = omp_get_wtime()
223 do isph = 1, params % nsph
224 call contract_grad_l(params, constants, isph, state % xs, &
225 & state % sgrid, workspace % tmp_vylm(:, 1), &
226 & workspace % tmp_vdylm(:, :, 1), workspace % tmp_vplm(:, 1), &
227 & workspace % tmp_vcos(:, 1), workspace % tmp_vsin(:, 1), &
229 call contract_grad_u(params, constants, isph, state % sgrid, &
230 & state % phi_grid, force(:, isph))
235 call zeta_grad(params, constants, state, e_cav, force)
237 finish_time = omp_get_wtime()
238 state % force_time = finish_time - start_time
256 & state, e_cav, force, dr, ddx_error)
258 type(ddx_params_type),
intent(in) :: params
259 type(ddx_constants_type),
intent(in) :: constants
260 type(ddx_workspace_type),
intent(inout) :: workspace
261 type(ddx_state_type),
intent(inout) :: state
262 type(ddx_error_type),
intent(inout) :: ddx_error
263 real(dp),
intent(in) :: e_cav(3, constants % ncav)
264 real(dp),
intent(inout) :: force(3, params % nsph)
265 real(dp),
intent(inout) :: dr(params % nsph)
267 real(dp) :: start_time, finish_time
271 if (ddx_error % flag .eq. 0)
continue
273 start_time = omp_get_wtime()
276 do isph = 1, params % nsph
277 call contract_grad_l(params, constants, isph, state % xs, &
278 & state % sgrid, workspace % tmp_vylm(:, 1), &
279 & workspace % tmp_vdylm(:, :, 1), workspace % tmp_vplm(:, 1), &
280 & workspace % tmp_vcos(:, 1), workspace % tmp_vsin(:, 1), &
281 & force(:, isph), dr(isph))
282 call contract_grad_u(params, constants, isph, state % sgrid, &
283 & state % phi_grid, force(:, isph), dr(isph))
288 call zeta_grad(params, constants, state, e_cav, force)
290 finish_time = omp_get_wtime()
291 state % force_time = finish_time - start_time
313 type(ddx_params_type),
intent(in) :: params
314 type(ddx_constants_type),
intent(in) :: constants
315 type(ddx_workspace_type),
intent(inout) :: workspace
316 type(ddx_state_type),
intent(inout) :: state
318 real(dp),
external :: ddot
319 integer :: icav, isph, igrid
322 if (
allocated(workspace % tmp_pot))
continue
325 call ddeval_grid_work(constants % nbasis, params % ngrid, params % nsph, &
326 & constants % vgrid, constants % vgrid_nbasis, one, state % s, zero, &
329 call ddcav_to_grid_work(params % ngrid, params % nsph, constants % ncav, &
330 & constants % icav_ia, constants % icav_ja, state % phi_cav, &
336 do isph = 1, params % nsph
337 do igrid = 1, params % ngrid
338 if (constants % ui(igrid, isph) .ne. zero)
then
340 state % zeta(icav) = constants % wgrid(igrid) * &
341 & constants % ui(igrid, isph) * ddot(constants % nbasis, &
342 & constants % vgrid(1, igrid), 1, state % s(1, isph), 1)
subroutine cosmo_guess(params, constants, workspace, state, ddx_error)
Do a guess for the primal ddCOSMO linear system.
subroutine cosmo_solve(params, constants, workspace, state, tol, ddx_error)
Solve the primal ddCOSMO linear system.
subroutine cosmo_solve_adjoint(params, constants, workspace, state, tol, ddx_error)
Solve the adjoint ddCOSMO linear system.
subroutine cosmo_solvation_force_dr_terms(params, constants, workspace, state, e_cav, force, dr, ddx_error)
Compute the solvation term of the forces and the radial derivatives (solute aspecific)....
subroutine cosmo_solvation_force_terms(params, constants, workspace, state, e_cav, force, ddx_error)
Compute the solvation term of the forces (solute aspecific). This must be summed to the solute specif...
subroutine cosmo_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.
subroutine cosmo_guess_adjoint(params, constants, workspace, state, ddx_error)
Do a guess for the adjoint ddCOSMO linear system.
subroutine cosmo_energy(constants, state, esolv, ddx_error)
Compute the ddCOSMO energy.
High-level subroutines for ddcosmo.
subroutine cosmo_derivative_setup(params, constants, workspace, state)
This routines precomputes the intermediates to be used in the evaluation of ddCOSMO analytical deriva...
Routines to build rhs (phi and psi)
Operators shared among ddX methods.
subroutine lx(params, constants, workspace, x, y, ddx_error)
Single layer operator matvec.
subroutine ldm1x(params, constants, workspace, x, y, ddx_error)
Diagonal preconditioning for Lx operator.
subroutine lstarx(params, constants, workspace, x, y, ddx_error)
Adjoint single layer operator matvec.