File: bxsf_write.c

package info (click to toggle)
c2x 2.42%2Bds-1
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 1,316 kB
  • sloc: ansic: 28,166; makefile: 31; sh: 1
file content (171 lines) | stat: -rw-r--r-- 5,153 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
/* Write an XCrysDen .bxsf Fermi surface file */

/* Copyright (c) 2020, 2025 MJ Rutter
 *
 * This program is free software; you can redistribute it and/or
 * modify it under the terms of the GNU General Public License as
 * published by the Free Software Foundation, either version 3
 * of the Licence, or (at your option) any later version.
 *
 * This program is distributed in the hope that it will be useful,
 * but WITHOUT ANY WARRANTY; without even the implied warranty of
 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
 * GNU General Public License for more details.
 *
 * You should have received a copy of the GNU General Public License
 * along with this program; if not, see http://www.gnu.org/licenses/
 */

#include<stdio.h>
#include<stdlib.h>
#include<math.h>

#include "c2xsf.h"

void bxsf_write(FILE* outfile, struct unit_cell *c, struct contents *m,
                struct es *elect, struct kpts *kp,
                struct symmetry *rs){
  int i,j,k,ii,jj,kk,ic,off,ind,nb,ns;
  double scale,lscale,frac1[3],chg,efermi;
  int grid[3],gpts,nbands,nspins;
  int *mapping;

  frac1[0]=frac1[1]=frac1[2]=0;
  scale=1;
  if (flags&AU) scale=1/H_eV;
  lscale=1;
  if (flags&AU) lscale=BOHR;
  
  if (!elect->eval)
    error_exit("No evalues found, so cannot write bxsf");
  
  /* Force kpoint co-ords to 0<=x<1 */

  for(i=0;i<kp->n;i++)
    for(j=0;j<3;j++)
      kp->kpts[i].frac[j]=fmod(kp->kpts[i].frac[j]+1.0,1.0);

  addabs(kp->kpts,kp->n,c->basis);

  if (!elect->e_fermi){
    efermi=NAN;
    chg=elect->nel;
    if (chg==0){
      for(i=0;i<m->n;i++)
	chg+=m->atoms[i].chg;
      if ((chg)&&(elect->charge)) chg-=*elect->charge;
      if (chg)
	fprintf(stderr,"Assuming nel=%lf from total ionic charge\n",chg);
    }
    if (chg) efermi=calc_efermi(elect,kp,chg);
    if (isnan(efermi)){
      efermi=0;
      fprintf(stderr,"Warning: no Fermi energy found\n");
    }
    else
      fprintf(stderr,"Warning: using estimate of Fermi energy\n");
  }
  else efermi=*elect->e_fermi;
  

  if (kp->mp){
    for(i=0;i<3;i++) grid[i]=kp->mp->grid[i];
  }
  else{
    grid_detect(kp,grid,frac1);
  }
 
  /* Check that the kpoint grid includes the gamma point */

  if (!aeq(vmod2(frac1),0.0))
      error_exit("Gamma point kpt must be included for .bxsf output");
  
  gpts=grid[0]*grid[1]*grid[2];

  if (gpts==0)
    error_exit("Unable to deduce k-point mesh parameters");

  mapping=malloc(gpts*sizeof(int));
  if (!mapping) error_exit("Malloc error for mapping");

  if (kgrid_expand(c,m,kp,rs,grid,frac1,mapping))
    error_exit("Failed to unfold kpoints to regular grid");
  
  
  fprintf(stderr,"nbspins=%d\n",elect->nbspins);
  /* Count bands */
  nbands=0;
  for(i=0;i<elect->nbands;i++)
    if (inrange(i+1,elect->band_range)) nbands++;
  /* Count spins */
  nspins=0;
  for(i=0;i<elect->nbspins;i++)
    if (inrange(i,elect->spin_range)) nspins++;
  
  fprintf(outfile,"BEGIN_INFO\n");
  fprintf(outfile,"# Generated by c2x\n");
  if (m->title) fprintf(outfile,"# %s\n",m->title);
  if (flags&AU)
    fprintf(outfile,"# Units: Hartrees and Bohr^-1\n");
  else
    fprintf(outfile,"# Units: eV and Angstrom^-1\n");
  if ((elect->nbspins>1)||(nbands!=elect->nbands)){
    fprintf(outfile,"#  Band number mapping\n");
    fprintf(outfile,"# number  label   band spin\n");
    i=1;
    for(ns=0;ns<elect->nbspins;ns++){
      if (!inrange(ns,elect->spin_range)) continue;
      for(nb=0;nb<elect->nbands;nb++){
        if (!inrange(nb+1,elect->band_range)) continue;
        fprintf(outfile,"# %3d     %4d   %3d    %1d\n",
                i,nb+1+elect->nbands*ns,nb+1,ns);
        i++;
      }
    }
  }
  fprintf(outfile,"  Fermi Energy: %.6f\n",efermi*scale);
  fprintf(outfile,"END_INFO\n\n");

  fprintf(outfile,"BEGIN_BLOCK_BANDGRID_3D\n");
  fprintf(outfile,"Comment\n");
  if (flags&AU)
    fprintf(outfile,"BEGIN_BANDGRID_3D_Ha\n");
  else
    fprintf(outfile,"BEGIN_BANDGRID_3D_eV\n");
  fprintf(outfile,"  %d\n",nbands*nspins);
  fprintf(outfile,"  %4d %4d %4d\n",grid[0]+1,grid[1]+1,grid[2]+1);
  fprintf(outfile,"  0.0 0.0 0.0\n");
  for(i=0;i<3;i++)
    fprintf(outfile,"  % 8f  % 8f  %8f\n",c->recip[i][0]*lscale,
            c->recip[i][1]*lscale,c->recip[i][2]*lscale);
  for(ns=0;ns<elect->nbspins;ns++){
    if (!inrange(ns,elect->spin_range)) continue;
    for(nb=0;nb<elect->nbands;nb++){
      if (!inrange(nb+1,elect->band_range)) continue;
      fprintf(outfile,"  BAND: %3d",nb+1+elect->nbands*ns);

 /* XCrysDen likes bizarre "general" grids */
      ic=0;
      for(i=0;i<=grid[0];i++){
        ii=i%grid[0];
        for(j=0;j<=grid[1];j++){
          jj=j%grid[1];
          off=jj*grid[2]+ii*grid[2]*grid[1];
          for(k=0;k<=grid[2];k++){
            kk=k%grid[2];
            ind=kk+off;
            if (ic%5==0) fprintf(outfile,"\n  ");
            fprintf(outfile,"  % 8f",elect->eval[nb+mapping[ind]*elect->nbands*elect->nbspins+elect->nbands*ns]*scale);
            ic++;
          }
        }
      }
      fprintf(outfile,"\n");
    }
  }
  fprintf(outfile,"END_BANDGRID_3D\n");
  fprintf(outfile,"END_BLOCK_BANDGRID_3D\n");

  free(mapping);
  print_bandwidths(elect,kp);
}