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
87subroutine cosmo_guess(params, constants, workspace, state, ddx_error)
89 type(ddx_params_type),
intent(in) :: params
90 type(ddx_constants_type),
intent(in) :: constants
91 type(ddx_workspace_type),
intent(inout) :: workspace
92 type(ddx_state_type),
intent(inout) :: state
93 type(ddx_error_type),
intent(inout) :: ddx_error
96 call ldm1x(params, constants, workspace, state % phi, state % xs, ddx_error)
111 type(ddx_params_type),
intent(in) :: params
112 type(ddx_constants_type),
intent(in) :: constants
113 type(ddx_workspace_type),
intent(inout) :: workspace
114 type(ddx_state_type),
intent(inout) :: state
115 type(ddx_error_type),
intent(inout) :: ddx_error
118 call ldm1x(params, constants, workspace, state % psi, state % s, ddx_error)
132subroutine cosmo_solve(params, constants, workspace, state, tol, ddx_error)
134 type(ddx_params_type),
intent(in) :: params
135 type(ddx_constants_type),
intent(in) :: constants
136 type(ddx_workspace_type),
intent(inout) :: workspace
137 type(ddx_state_type),
intent(inout) :: state
138 real(dp),
intent(in) :: tol
139 type(ddx_error_type),
intent(inout) :: ddx_error
141 real(dp) :: start_time, finish_time
143 state % xs_niter = params % maxiter
144 start_time = omp_get_wtime()
145 call jacobi_diis(params, constants, workspace, tol, state % phi, &
146 & state % xs, state % xs_niter, state % xs_rel_diff,
lx,
ldm1x, &
148 finish_time = omp_get_wtime()
149 state % xs_time = finish_time - start_time
166 type(ddx_params_type),
intent(in) :: params
167 type(ddx_constants_type),
intent(in) :: constants
168 type(ddx_workspace_type),
intent(inout) :: workspace
169 type(ddx_state_type),
intent(inout) :: state
170 real(dp),
intent(in) :: tol
171 type(ddx_error_type),
intent(inout) :: ddx_error
173 real(dp) :: start_time, finish_time
175 state % s_niter = params % maxiter
176 start_time = omp_get_wtime()
177 call jacobi_diis(params, constants, workspace, tol, state % psi, &
178 & state % s, state % s_niter, state % s_rel_diff,
lstarx,
ldm1x, &
180 finish_time = omp_get_wtime()
181 state % s_time = finish_time - start_time
183 state % q = state % s
202 & state, e_cav, force, ddx_error)
204 type(ddx_params_type),
intent(in) :: params
205 type(ddx_constants_type),
intent(in) :: constants
206 type(ddx_workspace_type),
intent(inout) :: workspace
207 type(ddx_state_type),
intent(inout) :: state
208 type(ddx_error_type),
intent(inout) :: ddx_error
209 real(dp),
intent(in) :: e_cav(3, constants % ncav)
210 real(dp),
intent(inout) :: force(3, params % nsph)
212 real(dp) :: start_time, finish_time
216 if (ddx_error % flag .eq. 0)
continue
218 start_time = omp_get_wtime()
221 do isph = 1, params % nsph
222 call contract_grad_l(params, constants, isph, state % xs, &
223 & state % sgrid, workspace % tmp_vylm(:, 1), &
224 & workspace % tmp_vdylm(:, :, 1), workspace % tmp_vplm(:, 1), &
225 & workspace % tmp_vcos(:, 1), workspace % tmp_vsin(:, 1), &
227 call contract_grad_u(params, constants, isph, state % sgrid, &
228 & state % phi_grid, force(:, isph))
233 call zeta_grad(params, constants, state, e_cav, force)
235 finish_time = omp_get_wtime()
236 state % force_time = finish_time - start_time
249 type(ddx_params_type),
intent(in) :: params
250 type(ddx_constants_type),
intent(in) :: constants
251 type(ddx_workspace_type),
intent(inout) :: workspace
252 type(ddx_state_type),
intent(inout) :: state
254 real(dp),
external :: ddot
255 integer :: icav, isph, igrid
258 if (
allocated(workspace % tmp_pot))
continue
261 call ddeval_grid_work(constants % nbasis, params % ngrid, params % nsph, &
262 & constants % vgrid, constants % vgrid_nbasis, one, state % s, zero, &
265 call ddcav_to_grid_work(params % ngrid, params % nsph, constants % ncav, &
266 & constants % icav_ia, constants % icav_ja, state % phi_cav, &
272 do isph = 1, params % nsph
273 do igrid = 1, params % ngrid
274 if (constants % ui(igrid, isph) .ne. zero)
then
276 state % zeta(icav) = constants % wgrid(igrid) * &
277 & constants % ui(igrid, isph) * ddot(constants % nbasis, &
278 & 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_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.