Merge lp:~albertog/siesta/4.0-spin into lp:siesta/4.0

Proposed by Alberto Garcia
Status: Merged
Approved by: Nick Papior
Approved revision: 560
Merged at revision: 562
Proposed branch: lp:~albertog/siesta/4.0-spin
Merge into: lp:siesta/4.0
Diff against target: 502 lines (+195/-88)
15 files modified
Docs/compatibility.tex (+5/-0)
Src/Makefile (+7/-7)
Src/compute_dm.F (+1/-9)
Src/mixer.F (+1/-1)
Src/mulliken.F (+29/-17)
Src/new_dm.F (+1/-43)
Src/print_spin.F90 (+9/-7)
Src/read_options.F90 (+1/-0)
Src/scfconvergence_test.F (+19/-0)
Src/siesta_options.F90 (+2/-1)
Src/state_init.F (+2/-0)
Util/DensityMatrix/README (+4/-0)
Util/DensityMatrix/dm_noncol_sign_flip4.f90 (+108/-0)
Util/DensityMatrix/makefile (+5/-2)
version.info (+1/-1)
To merge this branch: bzr merge lp:~albertog/siesta/4.0-spin
Reviewer Review Type Date Requested Status
Nick Papior Approve
Review via email: mp+342997@code.launchpad.net

Commit message

Spin monitoring in scf cycle. Clarification of sign convention for DM

  * Spin monitoring during scf cycle

  If the fdf variable 'Spin.In.Scf' is set to 'true', the size and
  components of the (total) spin polarization will be printed at every
  scf step. This is analogous to the 'Mulliken.In.Scf' feature.

  The default for the Spin.In.SCF fdf flag is 'true' for calculations
  involving spin.

  A single line of spin information is printed, prefixed by ' spin-moment:'.

  Note that now, for each scf step, the spin and/or mulliken analyses are
  printed *after* the line showing the energies.

  The initial spin polarization for every geometry iteration is also
  printed in 'state_init', replacing the similar calls in 'new_dm'.

  * Use explicitly the same DM_12 sign-convention as in SOC versions

  In the non-collinear case, the code is actually using internally the
  same sign convention for the "up-down" components of H and the DM as
  the spin-orbit-capable (SOC) versions (>= 4.1): The building of the
  dense Hamiltonian has the same form and the off-diagonal components
  of Vxc are identical. Routines diag2g and diag2k are actually using
  the same convention, but with DM* as an intermediate step.

  By changing the interface to the outside world, namely the
  correspondence between spin angles and the DM in new_dm, and the
  calculation of the spin (actually moment) components in spnvec, the
  code can now interoperate (i.e., reuse DMs) with the SOC versions.

  In order to reuse DM information produced by previous versions of
  the program in the non-collinear case, the sign of D(4) of the old
  DM has to be changed. This can be achieved with the
  dm_noncol_sign_flip4 program in Util/DensityMatrix.

To post a comment you must log in.
Revision history for this message
Nick Papior (nickpapior) wrote :

LGTM with the small changes. I think we should just delete the lines that you have commented out (by choosing better alternatives). No need to have them there... :)

review: Approve
Revision history for this message
Alberto Garcia (albertog) wrote :

I agree to remove the old sections.

> On 11 Apr 2018, at 11:41, Nick Papior <email address hidden> wrote:
>
> LGTM

Revision history for this message
Nick Papior (nickpapior) wrote :

Super! Great!

Preview Diff

[H/L] Next/Prev Comment, [J/K] Next/Prev File, [N/P] Next/Prev Hunk
1=== modified file 'Docs/compatibility.tex'
2--- Docs/compatibility.tex 2018-04-09 13:45:49 +0000
3+++ Docs/compatibility.tex 2018-04-11 09:36:55 +0000
4@@ -27,6 +27,11 @@
5 internal tables used to compute two-center integrals. Please see
6 \opt{Compat.Matel.NRTAB} for details.
7
8+ In non-collinear-spin calculations, the sign of the fourth spin
9+ component of the density matrix is reversed. To re-use old \code{DM} files,
10+ they can be converted using the
11+ \code{Util/DensityMatrix/dm_noncol_sign_flip4} program.
12+
13 \item[\emph{any} --- 4.0-b2] The following compatibility issues should be remarked when
14 comparing with any later version of \siesta.
15
16
17=== modified file 'Src/Makefile'
18--- Src/Makefile 2018-04-09 14:08:42 +0000
19+++ Src/Makefile 2018-04-11 09:36:55 +0000
20@@ -689,11 +689,10 @@
21 mulliken.o: alloc.o atmfuncs.o parallel.o parallelsubs.o precision.o
22 mulliken.o: siesta_cml.o
23 naefs.o: atmfuncs.o mneighb.o new_matel.o precision.o
24-new_dm.o: alloc.o atomlist.o m_energies.o m_fdf_global.o m_iodm.o m_mpi_utils.o
25-new_dm.o: m_sparse.o m_spin.o m_steps.o m_ts_global_vars.o m_ts_iodm.o
26-new_dm.o: m_ts_options.o normalize_dm.o parallel.o parallelsubs.o parsing.o
27-new_dm.o: precision.o siesta_geom.o siesta_options.o sparse_matrices.o sys.o
28-new_dm.o: units.o
29+new_dm.o: alloc.o atomlist.o m_energies.o m_fdf_global.o m_iodm.o m_sparse.o
30+new_dm.o: m_spin.o m_steps.o m_ts_global_vars.o m_ts_iodm.o m_ts_options.o
31+new_dm.o: normalize_dm.o parallel.o parallelsubs.o parsing.o precision.o
32+new_dm.o: siesta_geom.o siesta_options.o sparse_matrices.o sys.o units.o
33 new_matel.o: alloc.o errorf.o interpolation.o matel_registry.o parallel.o
34 new_matel.o: precision.o radfft.o siesta_options.o spher_harm.o sys.o
35 nlefsm.o: alloc.o atmfuncs.o mneighb.o new_matel.o parallel.o parallelsubs.o
36@@ -768,8 +767,9 @@
37 save_density_matrix.o: m_steps.o m_ts_global_vars.o m_ts_iodm.o precision.o
38 save_density_matrix.o: siesta_options.o sparse_matrices.o
39 savepsi.o: alloc.o parallel.o parallelsubs.o precision.o
40-scfconvergence_test.o: m_convergence.o m_energies.o m_wallclock.o parallel.o
41-scfconvergence_test.o: siesta_cml.o siesta_options.o units.o write_subs.o
42+scfconvergence_test.o: atomlist.o m_convergence.o m_energies.o m_spin.o
43+scfconvergence_test.o: m_wallclock.o parallel.o siesta_cml.o siesta_geom.o
44+scfconvergence_test.o: siesta_options.o sparse_matrices.o units.o write_subs.o
45 schecomm.o: alloc.o
46 setatomnodes.o: alloc.o parallel.o precision.o spatial.o sys.o
47 setspatial.o: alloc.o parallel.o precision.o spatial.o
48
49=== modified file 'Src/compute_dm.F'
50--- Src/compute_dm.F 2016-01-25 16:00:16 +0000
51+++ Src/compute_dm.F 2018-04-11 09:36:55 +0000
52@@ -51,6 +51,7 @@
53
54 ! e1>e2 to signal that we do not want DOS weights
55 real(dp), parameter :: e1 = 1.0_dp, e2 = -1.0_dp
56+ real(dp) :: dummy_qspin(8)
57
58 !-------------------------------------------------------------- BEGIN
59 call timer( 'compute_dm', 1 )
60@@ -147,15 +148,6 @@
61 endif
62 #endif
63
64-! Print populations at each SCF step if requested before mixing ......
65- if (muldeb) then
66- if (ionode) write (6,"(/a)")
67- & 'siesta: Mulliken populations (DM_out)'
68- call mulliken( mullipop, nspin, na_u, no_u, maxnh,
69- & numh, listhptr, listh, S, Dscf, isa,
70- & lasto, iaorb, iphorb )
71- endif
72-
73 ! Write orbital indexes. JMS Dec.2009
74
75 if (IOnode .and. iscf==1) then
76
77=== modified file 'Src/mixer.F'
78--- Src/mixer.F 2016-01-25 16:00:16 +0000
79+++ Src/mixer.F 2018-04-11 09:36:55 +0000
80@@ -94,7 +94,7 @@
81 ! entirely correct. It should be moved to the top,
82 ! or done somewhere else.
83
84- if (muldeb) then
85+ if (muldeb .and. (.not. mixH)) then
86 if (IONode)
87 & write (6,"(/a)")
88 & 'siesta: Mulliken populations after mixing'
89
90=== modified file 'Src/mulliken.F'
91--- Src/mulliken.F 2016-01-25 16:00:16 +0000
92+++ Src/mulliken.F 2018-04-11 09:36:55 +0000
93@@ -533,20 +533,26 @@
94 c integer ns : Number of components in spin density matrix
95 c real*8 qs(ns) : Spin density matrix elements with the convention
96 c is=1 => Q11; is=2 => Q22; is=3 => Real(Q12);
97-c is=4 => Imag(Q12)
98+c is=4 => -Imag(Q12) (*Note
99 c ******* Output *****************************************************
100 c real*8 qt : Total charge
101-c real*8 st : Total spin
102-c real*8 sv(3) : Spin vector
103+c real*8 st : Total spin moment
104+c real*8 sv(3) : Spin moment vector
105 c ********************************************************************
106-
107+!
108+! Spin magnetic **moments** in bohr magnetons:
109+! m_x = 2 Re{Q12} = 2*qs(3)
110+! m_y = -2 Im{Q12} = 2*qs(4)
111+! m_z = Q11-Q22 = qs(1)-qs(2)
112+!
113+! The spin moment for a spin angular momentum of 1/2 (hbar) is (very close to) 1 bohr magneton.
114+!
115 use precision
116
117 implicit none
118- integer ns
119- real(dp) qs(ns), qt, st, sv(3)
120- real(dp) cosph, costh, sinph, sinth, tiny
121- parameter ( tiny = 1.e-12_dp )
122+ integer, intent(in) :: ns
123+ real(dp), intent(in) :: qs(ns)
124+ real(dp), intent(out) :: qt, st, sv(3)
125
126 if (ns .eq. 1) then
127 qt = qs(1)
128@@ -561,15 +567,21 @@
129 sv(2) = 0.0_dp
130 sv(3) = st
131 elseif (ns .eq. 4) then
132- qt = qs(1) + qs(2)
133- st = sqrt( (qs(1)-qs(2))**2 + 4._dp*(qs(3)**2+qs(4)**2) )
134- costh = ( qs(1) - qs(2) ) / ( st + tiny )
135- sinth = sqrt( 1.0_dp - costh**2 )
136- cosph = qs(3) / ( sqrt( qs(3)**2 + qs(4)**2 ) + tiny )
137- sinph = -qs(4) / ( sqrt( qs(3)**2 + qs(4)**2 ) + tiny )
138- sv(1) = st * sinth * cosph
139- sv(2) = st * sinth * sinph
140- sv(3) = st * costh
141+ qt = qs(1) + qs(2)
142+ sv(1) = 2.0_dp * qs(3)
143+ sv(2) = 2.0_dp * qs(4)
144+ sv(3) = qs(1) - qs(2)
145+ st = sqrt(sv(1)**2+sv(2)**2+sv(3)**2)
146+!
147+! ... roundabout way ...
148+! st = sqrt( (qs(1)-qs(2))**2 + 4._dp*(qs(3)**2+qs(4)**2) )
149+! costh = ( qs(1) - qs(2) ) / ( st + tiny )
150+! sinth = sqrt( 1.0_dp - costh**2 )
151+! cosph = qs(3) / ( sqrt( qs(3)**2 + qs(4)**2 ) + tiny )
152+! sinph = qs(4) / ( sqrt( qs(3)**2 + qs(4)**2 ) + tiny )
153+! sv(1) = st * sinth * cosph
154+! sv(2) = st * sinth * sinph
155+! sv(3) = st * costh
156 else
157 write(6,*) 'spnvec: ERROR: invalid argument ns =', ns
158 return
159
160=== modified file 'Src/new_dm.F'
161--- Src/new_dm.F 2016-01-25 16:00:16 +0000
162+++ Src/new_dm.F 2018-04-11 09:36:55 +0000
163@@ -842,7 +842,7 @@
164 Dscf(ind,1) = (qio + spio * costh) / 2
165 Dscf(ind,2) = (qio - spio * costh) / 2
166 Dscf(ind,3) = spio * sinth * cosph / 2
167- Dscf(ind,4) = - spio * sinth * sinph / 2
168+ Dscf(ind,4) = spio * sinth * sinph / 2
169 else
170 Dscf(ind,1) = (qio + spio) / 2
171 Dscf(ind,2) = (qio - spio) / 2
172@@ -905,48 +905,6 @@
173
174 endif
175
176- call print_initial_spin()
177-
178-
179- CONTAINS
180-
181- subroutine print_initial_spin()
182- use m_mpi_utils, only: Globalize_sum
183- use sparse_matrices, only: S
184-
185- real(dp) :: qspin(nspin)
186- integer :: io, j, ispin, ind
187-#ifdef MPI
188- real(dp) :: qtmp(nspin)
189-#endif
190-
191-! Print spin polarization
192- if (nspin .ge. 2) then
193- do ispin = 1,nspin
194- qspin(ispin) = 0.0_dp
195- do io = 1,no_l
196- do j = 1,numh(io)
197- ind = listhptr(io)+j
198- qspin(ispin) = qspin(ispin) + Dscf(ind,ispin) * S(ind)
199- enddo
200- enddo
201- enddo
202-
203-#ifdef MPI
204-! Global reduction of spin components
205- call globalize_sum(qspin(1:nspin),qtmp(1:nspin))
206- qspin(1:nspin) = qtmp(1:nspin)
207-#endif
208- if (nspin .eq. 2) then
209- if (IOnode) then
210- write(6,'(/,a,f12.6)')
211- . 'initdm: Initial spin polarization (Qup-Qdown) =',
212- . qspin(1) - qspin(2)
213- endif
214- endif
215- endif
216- end subroutine print_initial_spin
217-
218 end subroutine initdm
219
220 end module m_new_dm
221
222=== modified file 'Src/print_spin.F90'
223--- Src/print_spin.F90 2017-03-05 16:38:31 +0000
224+++ Src/print_spin.F90 2018-04-11 09:36:55 +0000
225@@ -45,9 +45,11 @@
226 if (nspin .eq. 2) then
227
228 if (IOnode) then
229- write(6,'(/,a,f12.6)') &
230- 'siesta: Total spin polarization (Qup-Qdown) =', &
231- qspin(1) - qspin(2)
232+ Svec(1) = 0.0_dp
233+ Svec(2) = 0.0_dp
234+ Svec(3) = qspin(1) - qspin(2)
235+ Stot = Svec(3)
236+ write(6,'(5x,a,f10.5,2f10.1,f10.5)') 'spin moment: S , {S} = ', Stot, Svec
237 endif
238 if (cml_p) call cmlAddProperty(xf=mainXML, &
239 value=qspin(1)-qspin(2), dictref='siesta:stot', &
240@@ -57,10 +59,10 @@
241
242 call spnvec( nspin, qspin, qaux, Stot, Svec )
243 if (IONode) then
244- write(6,'(a,3f12.6)') &
245- 'siesta: Total spin polarization (x,y,z) = ', &
246- qspin(3)*2, qspin(4)*2, qspin(1)-qspin(2)
247- write(6,'(a,4f12.6)') 'siesta: S , {S} = ', Stot, Svec
248+! write(6,'(a,3f12.6)') &
249+! 'siesta: Total spin polarization (x,y,z) = ', &
250+! qspin(3)*2, -qspin(4)*2, qspin(1)-qspin(2)
251+ write(6,'(5x,a,4f10.5)') 'spin moment: S , {S} = ', Stot, Svec
252 if (cml_p) then
253 call cmlAddProperty(xf=mainXML, value=Stot, &
254 dictref='siesta:stot', units='siestaUnits:spin')
255
256=== modified file 'Src/read_options.F90'
257--- Src/read_options.F90 2018-04-05 07:11:04 +0000
258+++ Src/read_options.F90 2018-04-11 09:36:55 +0000
259@@ -1467,6 +1467,7 @@
260 dm_normalization_tol = fdf_get( 'DM.NormalizationTolerance',1.0d-5)
261 normalize_dm_during_scf= fdf_get( 'DM.NormalizeDuringSCF',.true.)
262 muldeb = fdf_get( 'MullikenInSCF' , .false.)
263+ spndeb = fdf_get( 'SpinInSCF' , (nspin>1) )
264 rijmin = fdf_get( 'WarningMinimumAtomicDistance', &
265 1.0_dp, 'Bohr' )
266 bornz = fdf_get( 'BornCharge' , .false. )
267
268=== modified file 'Src/scfconvergence_test.F'
269--- Src/scfconvergence_test.F 2016-01-25 16:00:16 +0000
270+++ Src/scfconvergence_test.F 2018-04-11 09:36:55 +0000
271@@ -23,6 +23,11 @@
272 use m_convergence, only: converger_t, tolerance
273 use m_convergence, only: add_value, is_converged
274
275+ use m_spin, only: nspin
276+ use sparse_matrices, only: Dscf, S, numh, listhptr, listh, maxnh
277+ use siesta_geom, only: na_u, isa
278+ use atomlist, only: no_u, lasto, iaorb, iphorb
279+
280 implicit none
281
282 integer :: iscf
283@@ -32,6 +37,7 @@
284 type(converger_t), intent(inout) :: conv_harris, conv_freeE
285 logical, intent(out) :: converged
286
287+ real(dp) :: dummy_qspin(8)
288 !------------------------------------------------------------------- BEGIN
289
290 ! convergence test
291@@ -59,6 +65,19 @@
292 call wallclock("-------------- end of scf step")
293 endif
294
295+! Print spin polarization at each SCF step if requested before mixing
296+ if (spndeb) then
297+ call print_spin(dummy_qspin)
298+ endif
299+! Print populations at each SCF step if requested before mixing ......
300+ if (muldeb) then
301+ if (ionode) write (6,"(/a)")
302+ & 'siesta: Mulliken populations (DM_out)'
303+ call mulliken( mullipop, nspin, na_u, no_u, maxnh,
304+ & numh, listhptr, listh, S, Dscf, isa,
305+ & lasto, iaorb, iphorb )
306+ endif
307+
308 converged = .false.
309
310 if (require_energy_convergence) then
311
312=== modified file 'Src/siesta_options.F90'
313--- Src/siesta_options.F90 2018-04-05 07:11:04 +0000
314+++ Src/siesta_options.F90 2018-04-11 09:36:55 +0000
315@@ -97,7 +97,8 @@
316
317 logical :: atmonly ! Set up pseudoatom information only?
318 logical :: harrisfun ! Use Harris functional?
319- logical :: muldeb ! Write Mulliken polpulations at every SCF step?
320+ logical :: muldeb ! Write Mulliken populations at every SCF step?
321+ logical :: spndeb ! Write spin-polarization information at every SCF step?
322 logical :: require_energy_convergence ! free Energy conv. to finish SCF iteration?
323 logical :: require_harris_convergence ! to finish SCF iteration?
324 logical :: broyden_optim ! Use Broyden method to optimize geometry?
325
326=== modified file 'Src/state_init.F'
327--- Src/state_init.F 2018-02-26 15:23:37 +0000
328+++ Src/state_init.F 2018-04-11 09:36:55 +0000
329@@ -70,6 +70,7 @@
330 logical :: auxchanged ! Auxiliary supercell changed?
331 logical :: folding, folding1
332 logical :: foundxv ! dummy for call to ioxv
333+ real(dp) :: dummy_qspin(8)
334 external :: madelung, timer
335 real(dp), external :: volcel
336
337@@ -332,6 +333,7 @@
338 ! Initialize density matrix
339 ! The resizing of Dscf is done inside new_dm
340 call new_dm( auxchanged, numh, listhptr, listh, Dscf)
341+ if (nspin > 1) call print_spin(dummy_qspin)
342
343 ! Check for size of Pulay auxiliary matrices
344 call init_pulay_arrays()
345
346=== modified file 'Util/DensityMatrix/README'
347--- Util/DensityMatrix/README 2009-08-20 09:25:37 +0000
348+++ Util/DensityMatrix/README 2018-04-11 09:36:55 +0000
349@@ -2,3 +2,7 @@
350 filter the contribution of certain orbitals to the DM.
351
352 Also, experimental octave scripts to read DM.nc and DMHS.nc files.
353+
354+To allow reuse of non-collinear DM files produced with versions <
355+4.0.2, use the program dm_noncol_sign_flip4.
356+
357
358=== added file 'Util/DensityMatrix/dm_noncol_sign_flip4.f90'
359--- Util/DensityMatrix/dm_noncol_sign_flip4.f90 1970-01-01 00:00:00 +0000
360+++ Util/DensityMatrix/dm_noncol_sign_flip4.f90 2018-04-11 09:36:55 +0000
361@@ -0,0 +1,108 @@
362+!
363+! Copyright (C) 1996-2018 The SIESTA group
364+! This file is distributed under the terms of the
365+! GNU General Public License: see COPYING in the top directory
366+! or http://www.gnu.org/copyleft/gpl.txt.
367+! See Docs/Contributors.txt for a list of contributors.
368+!
369+program dm_noncol_sign_flip4
370+
371+!
372+! Changes the sign of the "4" component of a DM (useful to re-use DM files created
373+! with versions prior to 4.0.2 in the non-collinear case).
374+!
375+! Usage: dm_noncol_sign_flip4
376+!
377+! The input file must be named "DM" (a symbolic link would do)
378+! The output file is called "DMSIGNFLIP4"
379+!
380+implicit none
381+
382+integer, parameter :: dp = selected_real_kind(14,100)
383+
384+integer :: norbs ! Number of atomic orbitals
385+integer :: nnzs ! Total number of orbital interactions
386+integer :: nspin ! Number of spins
387+
388+integer :: iog, ispin, iostat, j
389+
390+integer, dimension(:), allocatable :: numd
391+integer, dimension(:), allocatable :: row_pointer
392+integer, dimension(:), allocatable :: column
393+real(dp), dimension(:,:), allocatable :: dm
394+
395+!-----------------------------------------------------
396+
397+open(unit=1,file="DM",form="unformatted",status="old",action="read", &
398+ position="rewind",iostat=iostat)
399+if (iostat /= 0) then
400+ print *, "File DM cannot be opened"
401+ STOP
402+endif
403+
404+read(1) norbs, nspin
405+print *, "Norbs, nspin: ", norbs, nspin
406+
407+if (nspin /= 4) then
408+ print *, "The DM is not from a non-collinear spin calculation"
409+ STOP
410+endif
411+
412+allocate(numd(1:norbs),row_pointer(1:norbs))
413+
414+read(1) (numd(iog),iog=1,norbs)
415+nnzs = sum(numd(1:norbs))
416+!
417+allocate(column(1:nnzs),dm(1:nnzs,1:nspin))
418+
419+ ! Compute global row pointer
420+ row_pointer(1) = 0
421+ do iog=2,norbs
422+ row_pointer(iog) = row_pointer(iog-1) + numd(iog-1)
423+ enddo
424+
425+!
426+! Column information
427+!
428+do iog = 1, norbs
429+ read(1) (column(row_pointer(iog)+j),j=1,numd(iog))
430+enddo
431+
432+!
433+! DM
434+!
435+ do ispin = 1, nspin
436+ do iog = 1, norbs
437+ read(1) (dm(row_pointer(iog)+j,ispin),j=1,numd(iog))
438+ enddo
439+ enddo
440+
441+ close(1)
442+ !
443+ ! Change sign
444+ !
445+ dm(:,4) = -dm(:,4)
446+!
447+
448+open(unit=1,file="DMSIGNFLIP4",form="unformatted",status="new",action="write", &
449+ position="rewind",iostat=iostat)
450+if (iostat /= 0) then
451+ print *, "A file DMSIGNFLIP4 exists. Please remove or change its name"
452+ STOP
453+endif
454+
455+write(1) norbs, nspin
456+write(1) (numd(iog),iog=1,norbs)
457+do iog = 1, norbs
458+ write(1) (column(row_pointer(iog)+j),j=1,numd(iog))
459+enddo
460+do ispin = 1, nspin
461+ do iog = 1, norbs
462+ write(1) (dm(row_pointer(iog)+j,ispin),j=1,numd(iog))
463+ enddo
464+enddo
465+
466+close(1)
467+
468+
469+end program dm_noncol_sign_flip4
470
471=== modified file 'Util/DensityMatrix/makefile'
472--- Util/DensityMatrix/makefile 2016-01-25 16:00:16 +0000
473+++ Util/DensityMatrix/makefile 2018-04-11 09:36:55 +0000
474@@ -21,7 +21,7 @@
475 .SUFFIXES: .f .F .o .a .f90 .F90
476 #
477 default: dm2cdf cdf2dm
478-all: dm2cdf cdf2dm dmfilter
479+all: dm2cdf cdf2dm dmfilter dm_noncol_sign_flip4
480 #
481 OBJDIR=Obj
482 include ../../$(OBJDIR)/arch.make
483@@ -46,8 +46,11 @@
484 dmfilter: dmfilter.o
485 $(FC) $(LDFLAGS) -o $@ dmfilter.o
486 #
487+dm_noncol_sign_flip4: dm_noncol_sign_flip4.o
488+ $(FC) $(LDFLAGS) -o $@ dm_noncol_sign_flip4.o
489+#
490 clean:
491- rm -f *.o dm2cdf cdf2dm dmfilter *.o *.*d
492+ rm -f *.o dm2cdf cdf2dm dmfilter dm_noncol_sign_flip4 *.o *.*d
493 #
494
495
496
497=== modified file 'version.info'
498--- version.info 2018-04-09 14:08:42 +0000
499+++ version.info 2018-04-11 09:36:55 +0000
500@@ -1,1 +1,1 @@
501-siesta-4.0--561
502+siesta-4.0--561--spin-560

Subscribers

People subscribed via source and target branches