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 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394
|
/*****************************************************************************
Ver. 3.2 (April/01/2007)
OpenMX (Open source package for Material eXplorer) is a program package
for linear scaling density functional calculations of large-scale materials.
Almost time-consuming parts can be performed in O(N) operations where N is
the number of atoms (or basis orbitals). Thus, the program might be useful
for studies of large-scale materials.
The distribution of this program package follows the practice of
the GNU General Public Licence (GPL).
OpenMX is based on
* Local density and generalized gradient approximation (LDA, LSDA, GGA)
to the exchange-corellation term
* Norm-conserving pseudo potentials
* Variationally optimized pseudo atomic basis orbitals
* Solution of Poisson's equation using FFT
* Evaluation of two-center integrals using Fourier transformation
* Evaluation of three-center integrals using fixed real space grids
* Simple mixing, direct inversion in the interative subspace (DIIS),
and Guaranteed-reduction Pulay's methods for SCF calculations.
* Solution of the eigenvalue problem using O(N) methods
* ...
See also our website (http://staff.aist.go.jp/t-ozaki)
for recent developments.
**************************************************************
Copyright
Taisuke Ozaki
Present (Jan/01/2006) official address
RICS, National Institute of Advanced Industrial Science
and Technology (AIST), central 2, 1-1-1 Umezono, Tsukuba,
Ibaraki 305-8568, Japan
e-mail: t-ozaki@aist.go.jp
**************************************************************
*****************************************************************************/
/**********************************************************************
openmx.c:
openmx.c is the main routine of OpenMX.
Log of openmx.c:
5/Oct/2003 Released by T.Ozaki
***********************************************************************/
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
#include <string.h>
#include <time.h>
/* stat section */
#include <sys/types.h>
#include <sys/stat.h>
#include <unistd.h>
/* end stat section */
#include "openmx_common.h"
#ifdef nompi
#include "mimic_mpi.h"
#else
#include "mpi.h"
#endif
#ifdef TRAN
#include "tran_prototypes.h"
#include "tran_variables.h"
#endif
int main(int argc, char *argv[])
{
static int numprocs,myid;
static int MD_iter,i,j,po,ip;
static char fileE[YOUSO10] = ".ene";
static char fileDRC[YOUSO10] = ".md";
static char fileMemory[YOUSO10];
static char fileRestart[YOUSO10];
static char operate[200];
double TStime,TEtime;
static struct stat statbuf;
/* for idle CPUs */
int tag;
int complete;
MPI_Request request;
MPI_Status status;
/* MPI initialize */
mpi_comm_level1 = MPI_COMM_WORLD;
MPI_Init(&argc,&argv);
MPI_Comm_size(mpi_comm_level1,&numprocs);
MPI_Comm_rank(mpi_comm_level1,&myid);
Num_Procs = numprocs;
/* for measuring elapsed time */
dtime(&TStime);
/* check argv */
if (argc==1){
printf("\nCould not find an input file.\n\n");
MPI_Finalize();
exit(0);
}
/****************************************************
./openmx -maketest
making of *.out files in order to check whether
OpenMX normally runs on many platforms or not.
****************************************************/
if ( (argc==2 || argc==3) && strcmp(argv[1],"-maketest")==0){
Maketest(argc,argv);
exit(1);
}
/****************************************************
./openmx -runtest
check whether OpenMX normally runs on many platforms
or not by comparing the stored *.out and generated
*.out on your machine.
****************************************************/
else if ( strcmp(argv[1],"-runtest")==0){
Runtest("S",argc,argv);
}
/****************************************************
./openmx -maketestL
making of *.out files in order to check whether
OpenMX normally runs for relatively large systems
on many platforms or not
****************************************************/
if ( (argc==2 || argc==3) && strcmp(argv[1],"-maketestL")==0){
MaketestL(argc,argv);
exit(1);
}
/****************************************************
./openmx -runtestL
check whether OpenMX normally runs for relatively
large systems on many platforms or not by comparing
the stored *.out and generated *.out on your machine.
****************************************************/
else if (strcmp(argv[1],"-runtestL")==0){
Runtest("L",argc,argv);
}
/*******************************************************
check memory leak by monitoring actual used memory size
*******************************************************/
else if ( (argc==2 || argc==3) && strcmp(argv[1],"-mltest")==0){
Memory_Leak_test(argc,argv);
exit(1);
}
/*******************************************************
check consistency between analytic and numerical forces
*******************************************************/
else if ( (argc==3 || argc==4) && strcmp(argv[1],"-forcetest")==0){
if (strcmp(argv[2],"0")==0) force_flag = 0;
else if (strcmp(argv[2],"1")==0) force_flag = 1;
else if (strcmp(argv[2],"2")==0) force_flag = 2;
else if (strcmp(argv[2],"3")==0) force_flag = 3;
else if (strcmp(argv[2],"4")==0) force_flag = 4;
else if (strcmp(argv[2],"5")==0) force_flag = 5;
else if (strcmp(argv[2],"6")==0) force_flag = 6;
else if (strcmp(argv[2],"7")==0) force_flag = 7;
else if (strcmp(argv[2],"8")==0) force_flag = 8;
else {
printf("unsupported flag for -forcetest\n");
exit(1);
}
Force_test(argc,argv);
exit(1);
}
/* allocation of CompTime */
CompTime = (double**)malloc(sizeof(double*)*numprocs);
for (i=0; i<numprocs; i++){
CompTime[i] = (double*)malloc(sizeof(double)*20);
for (j=0; j<20; j++) CompTime[i][j] = 0.0;
}
if (myid==Host_ID){
printf("\n*******************************************************\n");
printf("*******************************************************\n");
printf(" Welcome to OpenMX Ver. %s \n",Version_OpenMX);
printf(" Copyright (C), 2002-2007, T.Ozaki \n");
printf(" OpenMX comes with ABSOLUTELY NO WARRANTY. \n");
printf(" This is free software, and you are welcome to \n");
printf(" redistribute it under the constitution of the GNU-GPL.\n");
printf("*******************************************************\n");
printf("*******************************************************\n\n");
}
Init_List_YOUSO();
remake_headfile = 0;
ScaleSize = 1.2;
/****************************************************
Read the input file
****************************************************/
/* setup CPU group */
setup_CPU_group(argv[1]);
if (myid>=atomnum) goto LAST_PROC; /* to cut off CPUs */
init_alloc_first();
CompTime[myid][1] = readfile(argv);
MPI_Barrier(mpi_comm_level1);
/* initialize PrintMemory routine */
sprintf(fileMemory,"%s%s.memory%i",filepath,filename,myid);
PrintMemory(fileMemory,0,"init");
PrintMemory_Fix();
/* initialize */
init();
fnjoint(filepath,filename,fileE);
fnjoint(filepath,filename,fileDRC);
/* check "-mltest2" mode */
po = 0;
if (myid==Host_ID){
for (i=0; i<argc; i++){
if ( strcmp(argv[i],"-mltest2")==0 ){
po = 1;
ip = i;
}
}
}
MPI_Bcast(&po, 1, MPI_INT, Host_ID, mpi_comm_level1);
MPI_Bcast(&ip, 1, MPI_INT, Host_ID, mpi_comm_level1);
if ( po==1 ) ML_flag = 1;
else ML_flag = 0;
/* check "-forcetest2" mode */
po = 0;
if (myid==Host_ID){
for (i=0; i<argc; i++){
if ( strcmp(argv[i],"-forcetest2")==0 ){
po = 1;
ip = i;
}
}
}
MPI_Bcast(&po, 1, MPI_INT, Host_ID, mpi_comm_level1);
MPI_Bcast(&ip, 1, MPI_INT, Host_ID, mpi_comm_level1);
if ( po==1 ){
force_flag = atoi(argv[ip+1]);
ForceConsistency_flag = 1;
}
/* check force consistency */
if (ForceConsistency_flag==1){
Check_Force(argv);
OutData(argv[1]);
Merge_LogFile(argv[1]);
Free_Arrays(0);
MPI_Finalize();
return 1;
}
/****************************************************
SCF-DFT calculations and MD and geometrical
optimization.
****************************************************/
MD_iter = 1;
do {
CompTime[myid][2] += truncation(MD_iter,Solver==6,1);
if (ML_flag==1 && myid==Host_ID) Get_VSZ(MD_iter);
#ifdef TRAN
if (Solver==4) {
TRAN_Calc_GridBound( mpi_comm_level1, atomnum, WhatSpecies, Spe_Atom_Cut1,
Ngrid1, Grid_Origin, Gxyz, tv, gtv, Right_tv );
/* output: TRAN_region[], TRAN_grid_bound */
}
#endif
CompTime[myid][3] += DFT(MD_iter,(MD_iter-1)%orbitalOpt_per_MDIter+1);
if (myid==Host_ID) iterout(MD_iter,MD_TimeStep*(MD_iter-1),fileE,fileDRC);
if (ML_flag==0) CompTime[myid][4] += MD_pac(MD_iter);
MD_iter++;
} while(MD_Opt_OK==0 && MD_iter<=MD_IterNumber);
#ifdef TRAN
if ( TRAN_output_hks ) {
/* left is dummy */
TRAN_RestartFile(mpi_comm_level1, "write","left",filepath,TRAN_hksoutfilename);
}
#endif
/****************************************************
calculate Voronoi charge
****************************************************/
if (Voronoi_Charge_flag==1) Voronoi_Charge();
/****************************************************
calculate Voronoi orbital magnetic moment
****************************************************/
if (Voronoi_OrbM_flag==1) Voronoi_Orbital_Moment();
/****************************************************
make an exchange splitting matrix
****************************************************/
if (ESM_flag==1) Make_ESM();
/****************************************************
Making of output files
****************************************************/
OutData(argv[1]);
/****************************************************
write connectivity, Hamiltonian, overlap, density
matrices, and etc. to a file, filename.scfout
****************************************************/
if (HS_fileout==1) SCF2File("write",argv[1]);
/* elapsed time */
dtime(&TEtime);
CompTime[myid][0] = TEtime - TStime;
Output_CompTime();
for (i=0; i<numprocs; i++){
free(CompTime[i]);
}
free(CompTime);
/* merge log files */
Merge_LogFile(argv[1]);
Free_Arrays(0);
LAST_PROC:
if (myid<atomnum) PrintMemory("total",0,"sum");
printf("\nThe calculation was normally finished. (proc=%3d)\n",myid);
MPI_Finalize();
return 0;
}
|