| [a0bcf1] | 1 | /** \file ions.c
|
|---|
| 2 | * Ionic Force, speed, coordinate and energy calculation.
|
|---|
| 3 | * Herein are all the routines for the molecular dynamics calculations such as
|
|---|
| 4 | * reading of configuration files IonsInitRead() and
|
|---|
| 5 | * free'ing RemoveIonsRead(),
|
|---|
| 6 | * summing up forces CalculateIonForce(),
|
|---|
| 7 | * correcting for temperature CorrectVelocities(),
|
|---|
| 8 | * or for fixation center of balance CorrectForces(),
|
|---|
| 9 | * outputting them to file OutputIonForce(),
|
|---|
| 10 | * calculating kinetic CalculateEnergyIonsU(), ewald CalculateEwald() energies,
|
|---|
| 11 | * moving ion UpdateIonsR() (position) UpdateIonsU() (speed),
|
|---|
| 12 | * scaling to match some temperature ScaleTemp()
|
|---|
| 13 | * are gathered.
|
|---|
| 14 | * Finally, needed for structure optimization GetOuterStop() evaluates the stop
|
|---|
| 15 | * condition of a minimal mean force acting on each ion.
|
|---|
| 16 | *
|
|---|
| 17 | Project: ParallelCarParrinello
|
|---|
| 18 | \author Jan Hamaekers
|
|---|
| 19 | \date 2000
|
|---|
| 20 |
|
|---|
| 21 | File: ions.c
|
|---|
| 22 | $Id: ions.c,v 1.34 2007/02/05 16:15:27 foo Exp $
|
|---|
| 23 | */
|
|---|
| 24 |
|
|---|
| [79290f] | 25 | #ifdef HAVE_CONFIG_H
|
|---|
| 26 | #include <config.h>
|
|---|
| 27 | #endif
|
|---|
| 28 |
|
|---|
| [a0bcf1] | 29 | #include <stdlib.h>
|
|---|
| 30 | #include <stdio.h>
|
|---|
| 31 | #include <stddef.h>
|
|---|
| 32 | #include <math.h>
|
|---|
| 33 | #include <string.h>
|
|---|
| 34 | #include <unistd.h>
|
|---|
| [961b75] | 35 | #include <gsl/gsl_randist.h>
|
|---|
| [a0bcf1] | 36 | #include "mymath.h"
|
|---|
| 37 |
|
|---|
| 38 | // use double precision fft when we have it
|
|---|
| 39 | #ifdef HAVE_DFFTW_H
|
|---|
| 40 | #include "dfftw.h"
|
|---|
| 41 | #else
|
|---|
| 42 | #include "fftw.h"
|
|---|
| 43 | #endif
|
|---|
| 44 |
|
|---|
| 45 | #include "data.h"
|
|---|
| 46 | #include "errors.h"
|
|---|
| 47 | #include "helpers.h"
|
|---|
| 48 | #include "mymath.h"
|
|---|
| 49 | #include "ions.h"
|
|---|
| 50 | #include "init.h"
|
|---|
| 51 |
|
|---|
| [b74e5e] | 52 | /** Readin of Thermostat related values from parameter file.
|
|---|
| 53 | * \param *P Problem at hand
|
|---|
| 54 | * \param *source parameter file
|
|---|
| 55 | */
|
|---|
| 56 | void InitThermostats(struct Problem *P, FILE *source)
|
|---|
| 57 | {
|
|---|
| 58 | struct FileData *F = &P->Files;
|
|---|
| 59 | struct Ions *I = &P->Ion;
|
|---|
| 60 | char *thermo = MallocString(12, "IonsInitRead: thermo");
|
|---|
| 61 | int j;
|
|---|
| 62 |
|
|---|
| 63 | // initialize the naming stuff and all
|
|---|
| 64 | I->ThermostatImplemented = (int *) Malloc((MaxThermostats)*(sizeof(int)), "IonsInitRead: *ThermostatImplemented");
|
|---|
| 65 | I->ThermostatNames = (char **) Malloc((MaxThermostats)*(sizeof(char *)), "IonsInitRead: *ThermostatNames");
|
|---|
| 66 | for (j=0;j<MaxThermostats;j++)
|
|---|
| 67 | I->ThermostatNames[j] = (char *) MallocString(12*(sizeof(char)), "IonsInitRead: ThermostatNames[]");
|
|---|
| 68 | strcpy(I->ThermostatNames[0],"None");
|
|---|
| 69 | I->ThermostatImplemented[0] = 1;
|
|---|
| 70 | strcpy(I->ThermostatNames[1],"Woodcock");
|
|---|
| 71 | I->ThermostatImplemented[1] = 1;
|
|---|
| 72 | strcpy(I->ThermostatNames[2],"Gaussian");
|
|---|
| 73 | I->ThermostatImplemented[2] = 1;
|
|---|
| 74 | strcpy(I->ThermostatNames[3],"Langevin");
|
|---|
| 75 | I->ThermostatImplemented[3] = 1;
|
|---|
| 76 | strcpy(I->ThermostatNames[4],"Berendsen");
|
|---|
| 77 | I->ThermostatImplemented[4] = 1;
|
|---|
| 78 | strcpy(I->ThermostatNames[5],"NoseHoover");
|
|---|
| 79 | I->ThermostatImplemented[5] = 1;
|
|---|
| 80 | I->Thermostat = -1;
|
|---|
| 81 | if (F->MeOutMes) fprintf(F->TemperatureFile, "# %s\t%s\t%s --- ","Step", "ActualTemp(1)", "ActualTemp(2)");
|
|---|
| 82 |
|
|---|
| 83 | // read desired Thermostat from file along with needed additional parameters
|
|---|
| 84 | if (ParseForParameter(P->Call.out[ReadOut],source,"Thermostat", 0, 1, 1, string_type, thermo, 1, optional)) {
|
|---|
| 85 | if (strcmp(thermo, I->ThermostatNames[0]) == 0) { // None
|
|---|
| 86 | if (I->ThermostatImplemented[0] == 1) {
|
|---|
| 87 | I->Thermostat = None;
|
|---|
| 88 | } else {
|
|---|
| 89 | if (P->Par.me == 0)
|
|---|
| 90 | fprintf(stderr, "(%i) Warning: %s thermostat not implemented, falling back to None.", P->Par.me, I->ThermostatNames[0]);
|
|---|
| 91 | I->Thermostat = None;
|
|---|
| 92 | }
|
|---|
| 93 | } else if (strcmp(thermo, I->ThermostatNames[1]) == 0) { // Woodcock
|
|---|
| 94 | if (I->ThermostatImplemented[1] == 1) {
|
|---|
| 95 | I->Thermostat = Woodcock;
|
|---|
| 96 | ParseForParameter(P->Call.out[ReadOut],source,"Thermostat", 0, 2, 1, int_type, &P->R.ScaleTempStep, 1, critical); // read scaling frequency
|
|---|
| 97 | } else {
|
|---|
| 98 | if (P->Par.me == 0)
|
|---|
| 99 | fprintf(stderr, "(%i) Warning: %s thermostat not implemented, falling back to None.", P->Par.me, I->ThermostatNames[1]);
|
|---|
| 100 | I->Thermostat = None;
|
|---|
| 101 | }
|
|---|
| 102 | } else if (strcmp(thermo, I->ThermostatNames[2]) == 0) { // Gaussian
|
|---|
| 103 | if (I->ThermostatImplemented[2] == 1) {
|
|---|
| 104 | I->Thermostat = Gaussian;
|
|---|
| 105 | ParseForParameter(P->Call.out[ReadOut],source,"Thermostat", 0, 2, 1, int_type, &P->R.ScaleTempStep, 1, critical); // read collision rate
|
|---|
| 106 | } else {
|
|---|
| 107 | if (P->Par.me == 0)
|
|---|
| 108 | fprintf(stderr, "(%i) Warning: %s thermostat not implemented, falling back to None.", P->Par.me, I->ThermostatNames[2]);
|
|---|
| 109 | I->Thermostat = None;
|
|---|
| 110 | }
|
|---|
| 111 | } else if (strcmp(thermo, I->ThermostatNames[3]) == 0) { // Langevin
|
|---|
| 112 | if (I->ThermostatImplemented[3] == 1) {
|
|---|
| 113 | I->Thermostat = Langevin;
|
|---|
| 114 | ParseForParameter(P->Call.out[ReadOut],source,"Thermostat", 0, 2, 1, double_type, &P->R.TempFrequency, 1, critical); // read gamma
|
|---|
| 115 | if (ParseForParameter(P->Call.out[ReadOut],source,"Thermostat", 0, 3, 1, double_type, &P->R.alpha, 1, optional)) {
|
|---|
| 116 | if (P->Par.me == 0)
|
|---|
| 117 | fprintf(stderr,"(%i) Extended Stochastic Thermostat detected with interpolation coefficient %lg\n", P->Par.me, P->R.alpha);
|
|---|
| 118 | } else {
|
|---|
| 119 | P->R.alpha = 1.;
|
|---|
| 120 | }
|
|---|
| 121 | } else {
|
|---|
| 122 | if (P->Par.me == 0)
|
|---|
| 123 | fprintf(stderr, "(%i) Warning: %s thermostat not implemented, falling back to None.", P->Par.me, I->ThermostatNames[3]);
|
|---|
| 124 | I->Thermostat = None;
|
|---|
| 125 | }
|
|---|
| 126 | } else if (strcmp(thermo, I->ThermostatNames[4]) == 0) { // Berendsen
|
|---|
| 127 | if (I->ThermostatImplemented[4] == 1) {
|
|---|
| 128 | I->Thermostat = Berendsen;
|
|---|
| 129 | ParseForParameter(P->Call.out[ReadOut],source,"Thermostat", 0, 2, 1, double_type, &P->R.TempFrequency, 1, critical); // read \tau_T
|
|---|
| 130 | } else {
|
|---|
| 131 | if (P->Par.me == 0)
|
|---|
| 132 | fprintf(stderr, "(%i) Warning: %s thermostat not implemented, falling back to None.", P->Par.me, I->ThermostatNames[4]);
|
|---|
| 133 | I->Thermostat = None;
|
|---|
| 134 | }
|
|---|
| 135 | } else if (strcmp(thermo, I->ThermostatNames[5]) == 0) { // Nose-Hoover
|
|---|
| 136 | if (I->ThermostatImplemented[5] == 1) {
|
|---|
| 137 | I->Thermostat = NoseHoover;
|
|---|
| 138 | ParseForParameter(P->Call.out[ReadOut],source,"Thermostat", 0, 2, 1, double_type, &P->R.HooverMass, 1, critical); // read Hoovermass
|
|---|
| 139 | P->R.alpha = 0.;
|
|---|
| 140 | } else {
|
|---|
| 141 | if (P->Par.me == 0)
|
|---|
| 142 | fprintf(stderr, "(%i) Warning: %s thermostat not implemented, falling back to None.", P->Par.me, I->ThermostatNames[5]);
|
|---|
| 143 | I->Thermostat = None;
|
|---|
| 144 | }
|
|---|
| [79290f] | 145 | } else {
|
|---|
| [b74e5e] | 146 | if (P->Par.me == 0)
|
|---|
| 147 | fprintf(stderr,"(%i) Warning: thermostat name was not understood!\n", P->Par.me);
|
|---|
| 148 | I->Thermostat = None;
|
|---|
| 149 | }
|
|---|
| 150 | } else {
|
|---|
| 151 | if ((P->R.MaxOuterStep > 0) && (P->R.TargetTemp != 0))
|
|---|
| 152 | debug(P, "No thermostat chosen despite finite temperature MD, falling back to None.");
|
|---|
| 153 | I->Thermostat = None;
|
|---|
| 154 | }
|
|---|
| 155 | if (F->MeOutMes) fprintf(F->TemperatureFile, "%s Thermostat\n",I->ThermostatNames[(int)I->Thermostat]);
|
|---|
| 156 | Free(thermo, "InitThermostats: thermo");
|
|---|
| 157 | }
|
|---|
| 158 |
|
|---|
| [961b75] | 159 | /** Parses for a given keyword the ion coordinates and velocites.
|
|---|
| 160 | * \param *P Problem at hand, contains pointer to Ion and IonType structure
|
|---|
| 161 | * \param *source
|
|---|
| 162 | * \param *keyword keyword string
|
|---|
| 163 | * \param repetition which repeated version of the keyword shall be read by ReadParameters()
|
|---|
| 164 | * \param relative whether atom coordinates are relative to unit cell size or absolute
|
|---|
| 165 | * \param Bohr conversion factor to atomiclength (e.g. 1.0 if coordinates in bohr radius, 0.529... if given in Angstrom)
|
|---|
| 166 | * \param sigma variance of maxwell-boltzmann distribution
|
|---|
| 167 | * \param *R where to put the read coordinates
|
|---|
| 168 | * \param *U where to put the read velocities
|
|---|
| 169 | * \param *IMT whether its MoveType#FixedIon or not
|
|---|
| 170 | * \return if 1 - success, 0 - failure
|
|---|
| 171 | */
|
|---|
| 172 | int ParseIonsCoordinatesAndVelocities(struct Problem *P, FILE *source, char *keyword, int repetition, int relative, double sigma, double *R, double *U, int *IMT)
|
|---|
| 173 | {
|
|---|
| 174 | //struct RunStruct *Run = &P->R;
|
|---|
| 175 | struct Ions *I = &P->Ion;
|
|---|
| 176 | struct Lattice *L = &P->Lat;
|
|---|
| 177 | int k;
|
|---|
| 178 | double R_trafo[3];
|
|---|
| 179 | gsl_rng * r;
|
|---|
| 180 | const gsl_rng_type * T;
|
|---|
| 181 | double Bohr = (I->IsAngstroem) ? ANGSTROEMINBOHRRADIUS : 1.;
|
|---|
| 182 |
|
|---|
| 183 | gsl_rng_env_setup();
|
|---|
| 184 | T = gsl_rng_default;
|
|---|
| 185 | r = gsl_rng_alloc (T);
|
|---|
| 186 | if (ParseForParameter(P->Call.out[ReadOut],source, keyword, 0, 1, 1, double_type, &R_trafo[0], repetition, optional) &&
|
|---|
| 187 | ParseForParameter(P->Call.out[ReadOut],source, keyword, 0, 2, 1, double_type, &R_trafo[1], repetition, optional) &&
|
|---|
| 188 | ParseForParameter(P->Call.out[ReadOut],source, keyword, 0, 3, 1, double_type, &R_trafo[2], repetition, optional) &&
|
|---|
| 189 | ParseForParameter(P->Call.out[ReadOut],source, keyword, 0, 4, 1, int_type, IMT, repetition, optional)) {
|
|---|
| 190 | // if read successful, then try also parsing velocities
|
|---|
| 191 | if ((ParseForParameter(P->Call.out[ReadOut],source, keyword, 0, 5, 1, double_type, &U[0], repetition, optional)) &&
|
|---|
| 192 | (ParseForParameter(P->Call.out[ReadOut],source, keyword, 0, 6, 1, double_type, &U[1], repetition, optional)) &&
|
|---|
| 193 | (ParseForParameter(P->Call.out[ReadOut],source, keyword, 0, 7, 1, double_type, &U[2], repetition, optional))) {
|
|---|
| 194 | //nothing
|
|---|
| 195 | } else { // if no velocities specified, set to maxwell-boltzmann distributed
|
|---|
| 196 | if ((I->Thermostat != None)) {
|
|---|
| 197 | if (P->Par.me == 0) fprintf(stderr, "(%i) Missing velocities of ion %s: Maxwell-Boltzmann with %lg\n", P->Par.me, keyword, sigma);
|
|---|
| 198 | for (k=0;k<NDIM;k++)
|
|---|
| 199 | U[k] = gsl_ran_gaussian (r, sigma);
|
|---|
| 200 | } else {
|
|---|
| 201 | for (k=0;k<NDIM;k++)
|
|---|
| 202 | U[k] = 0.;
|
|---|
| 203 | }
|
|---|
| 204 | }
|
|---|
| 205 | // change if coordinates were relative
|
|---|
| 206 | if (relative) {
|
|---|
| 207 | //fprintf(stderr,"(%i)Ion coordinates are relative %i ... \n",P->Par.me, relative);
|
|---|
| 208 | RMat33Vec3(R, L->RealBasis, R_trafo); // multiply with real basis
|
|---|
| 209 | } else {
|
|---|
| 210 | for (k=0;k<NDIM;k++) // otherweise just copy to vector
|
|---|
| 211 | R[k] = R_trafo[k];
|
|---|
| 212 | }
|
|---|
| 213 | // check if given MoveType is valid
|
|---|
| 214 | if (*IMT < 0 || *IMT >= MaxIonMoveType) {
|
|---|
| 215 | fprintf(stderr,"Bad Ion MoveType set to MoveIon for %s\n", keyword);
|
|---|
| 216 | *IMT = 0;
|
|---|
| 217 | }
|
|---|
| 218 | if (!relative) {
|
|---|
| 219 | //fprintf(stderr,"converting with %lg\n", Bohr);
|
|---|
| 220 | SM(R, Bohr, NDIM); // convert angstroem to bohrradius
|
|---|
| 221 | SM(U, Bohr, NDIM); // convert angstroem to bohrradius
|
|---|
| 222 | }
|
|---|
| 223 | // finally periodicly correct coordinates to remain in unit cell
|
|---|
| 224 | RMat33Vec3(R_trafo, L->InvBasis, R);
|
|---|
| 225 | for (k=0; k <NDIM; k++) { // periodic correction until within unit cell
|
|---|
| 226 | while (R_trafo[k] < 0)
|
|---|
| 227 | R_trafo[k] += 1.0;
|
|---|
| 228 | while (R_trafo[k] >= 1.0)
|
|---|
| 229 | R_trafo[k] -= 1.0;
|
|---|
| 230 | }
|
|---|
| 231 | RMat33Vec3(R, L->RealBasis, R_trafo);
|
|---|
| 232 | // print coordinates if desired
|
|---|
| 233 | if ((P->Call.out[ReadOut])) {
|
|---|
| 234 | fprintf(stderr,"(%i) coordinates of Ion %s: (x,y,z) = (",P->Par.me, keyword);
|
|---|
| 235 | for(k=0;k<NDIM;k++) {
|
|---|
| 236 | fprintf(stderr,"%lg ",R[k]);
|
|---|
| 237 | if (k != NDIM-1) fprintf(stderr,", ");
|
|---|
| 238 | else fprintf(stderr,")\n");
|
|---|
| 239 | }
|
|---|
| 240 | }
|
|---|
| 241 | gsl_rng_free (r);
|
|---|
| 242 | return 1;
|
|---|
| 243 | } else {
|
|---|
| 244 | for (k=0;k<NDIM;k++) {
|
|---|
| 245 | U[k] = R[k] = 0.;
|
|---|
| 246 | }
|
|---|
| 247 | *IMT = 0;
|
|---|
| 248 | gsl_rng_free (r);
|
|---|
| 249 | return 0;
|
|---|
| 250 | }
|
|---|
| 251 | }
|
|---|
| 252 |
|
|---|
| [a0bcf1] | 253 | /** Readin of the ion section of parameter file.
|
|---|
| 254 | * Among others the following paramaters are read from file:
|
|---|
| 255 | * Ions::MaxTypes, IonType::Max_IonsOfType, IonType::Z, IonType::R,
|
|---|
| 256 | * IonType:IMT, ...
|
|---|
| 257 | * \param P Problem at hand (containing Ions and Lattice structures)
|
|---|
| 258 | * \param *source file pointer
|
|---|
| 259 | */
|
|---|
| 260 | void IonsInitRead(struct Problem *P, FILE *source) {
|
|---|
| 261 | struct Ions *I = &P->Ion;
|
|---|
| 262 | //struct RunStruct *Run = &P->R;
|
|---|
| [961b75] | 263 | int i,it,j,IMT,d,me, count;
|
|---|
| [a0bcf1] | 264 | int relative; // whether read-in coordinate are relative (1) or absolute
|
|---|
| [460e95] | 265 | double Bohr = ANGSTROEMINBOHRRADIUS; /* Angstroem */
|
|---|
| 266 | double sigma, R[3], U[3];
|
|---|
| [961b75] | 267 | struct IonConstrained *ptr = NULL;
|
|---|
| 268 |
|
|---|
| [a0bcf1] | 269 | I->Max_TotalIons = 0;
|
|---|
| 270 | I->TotalZval = 0;
|
|---|
| 271 | MPI_Comm_rank(MPI_COMM_WORLD, &me);
|
|---|
| [961b75] | 272 | ParseForParameter(P->Call.out[ReadOut],source,"RCut", 0, 1, 1, double_type, &I->R_cut, 1, critical);
|
|---|
| 273 | ParseForParameter(P->Call.out[ReadOut],source,"IsAngstroem", 0, 1, 1, int_type, &I->IsAngstroem, 1, critical);
|
|---|
| [460e95] | 274 | if (!I->IsAngstroem) {
|
|---|
| [961b75] | 275 | Bohr = 1.0;
|
|---|
| 276 | } else { // adapt RealBasis (is in Angstroem as well)
|
|---|
| [460e95] | 277 | SM(P->Lat.RealBasis, Bohr, NDIM_NDIM);
|
|---|
| 278 | RMatReci3(P->Lat.InvBasis, P->Lat.RealBasis);
|
|---|
| 279 | }
|
|---|
| [961b75] | 280 | ParseForParameter(P->Call.out[ReadOut],source,"MaxTypes", 0, 1, 1, int_type, &I->Max_Types, 1, critical);
|
|---|
| [a0bcf1] | 281 | I->I = (struct IonType *) Malloc(sizeof(struct IonType)*I->Max_Types, "IonsInitRead: IonType");
|
|---|
| [961b75] | 282 | if (!ParseForParameter(P->Call.out[ReadOut],source,"RelativeCoord", 0, 1, 1, int_type, &relative, 1, optional))
|
|---|
| [a0bcf1] | 283 | relative = 0;
|
|---|
| [961b75] | 284 | if (!ParseForParameter(P->Call.out[ReadOut],source,"StructOpt", 0, 1, 1, int_type, &I->StructOpt, 1, optional))
|
|---|
| [a0bcf1] | 285 | I->StructOpt = 0;
|
|---|
| [961b75] | 286 |
|
|---|
| [b74e5e] | 287 | InitThermostats(P, source);
|
|---|
| [961b75] | 288 |
|
|---|
| [a0bcf1] | 289 | /* Ions Data */
|
|---|
| 290 | I->RLatticeVec = NULL;
|
|---|
| 291 | I->TotalMass = 0;
|
|---|
| 292 | char *free_name, *name;
|
|---|
| [c510a7] | 293 | name = free_name = Malloc(MAXSTRINGSIZE*sizeof(char),"IonsInitRead: Name");
|
|---|
| [a0bcf1] | 294 | for (i=0; i < I->Max_Types; i++) {
|
|---|
| 295 | sprintf(name,"Ion_Type%i",i+1);
|
|---|
| 296 | I->I[i].corecorr = NotCoreCorrected;
|
|---|
| [c510a7] | 297 | I->I[i].Name = MallocString(MAXSTRINGSIZE, "IonsInitRead: Name");
|
|---|
| 298 | I->I[i].Symbol = MallocString(MAXSTRINGSIZE, "IonsInitRead: Symbol");
|
|---|
| [961b75] | 299 | ParseForParameter(P->Call.out[ReadOut],source, name, 0, 1, 1, int_type, &I->I[i].Max_IonsOfType, 1, critical);
|
|---|
| 300 | ParseForParameter(P->Call.out[ReadOut],source, name, 0, 2, 1, int_type, &I->I[i].Z, 1, critical);
|
|---|
| 301 | ParseForParameter(P->Call.out[ReadOut],source, name, 0, 3, 1, double_type, &I->I[i].rgauss, 1, critical);
|
|---|
| 302 | ParseForParameter(P->Call.out[ReadOut],source, name, 0, 4, 1, int_type, &I->I[i].l_max, 1, critical);
|
|---|
| 303 | ParseForParameter(P->Call.out[ReadOut],source, name, 0, 5, 1, int_type, &I->I[i].l_loc, 1, critical);
|
|---|
| 304 | ParseForParameter(P->Call.out[ReadOut],source, name, 0, 6, 1, double_type, &I->I[i].IonMass, 1, critical);
|
|---|
| 305 | if(!ParseForParameter(P->Call.out[ReadOut],source, name, 0, 7, 1, string_type, &I->I[i].Name[0], 1, optional))
|
|---|
| [a0bcf1] | 306 | sprintf(&I->I[i].Name[0],"type%i",i);
|
|---|
| [961b75] | 307 | if(!ParseForParameter(P->Call.out[ReadOut],source, name, 0, 8, 1, string_type, &I->I[i].Symbol[0], 1, optional))
|
|---|
| [a0bcf1] | 308 | sprintf(&I->I[i].Symbol[0],"t%i",i);
|
|---|
| 309 | if (P->Call.out[ReadOut]) fprintf(stderr,"(%i) Element name: %s, symbol: %s\n",P->Par.me, I->I[i].Name, I->I[i].Symbol);
|
|---|
| [961b75] | 310 | I->I[i].moment = (double **) Malloc(sizeof(double *) * I->I[i].Max_IonsOfType, "IonsInitRead: moment");
|
|---|
| [a0bcf1] | 311 | I->I[i].sigma = (double **) Malloc(sizeof(double *) * I->I[i].Max_IonsOfType, "IonsInitRead: sigma");
|
|---|
| 312 | I->I[i].sigma_rezi = (double **) Malloc(sizeof(double *) * I->I[i].Max_IonsOfType, "IonsInitRead: sigma_rezi");
|
|---|
| 313 | I->I[i].sigma_PAS = (double **) Malloc(sizeof(double *) * I->I[i].Max_IonsOfType, "IonsInitRead: sigma_PAS");
|
|---|
| 314 | I->I[i].sigma_rezi_PAS = (double **) Malloc(sizeof(double *) * I->I[i].Max_IonsOfType, "IonsInitRead: sigma_rezi_PAS");
|
|---|
| 315 | for (it=0;it<I->I[i].Max_IonsOfType;it++) {
|
|---|
| [961b75] | 316 | I->I[i].moment[it] = (double *) Malloc(sizeof(double) * NDIM*NDIM, "IonsInitRead: moment");
|
|---|
| [a0bcf1] | 317 | I->I[i].sigma[it] = (double *) Malloc(sizeof(double) * NDIM*NDIM, "IonsInitRead: sigma");
|
|---|
| 318 | I->I[i].sigma_rezi[it] = (double *) Malloc(sizeof(double) * NDIM*NDIM, "IonsInitRead: sigma_rezi");
|
|---|
| 319 | I->I[i].sigma_PAS[it] = (double *) Malloc(sizeof(double) * NDIM, "IonsInitRead: sigma_PAS");
|
|---|
| 320 | I->I[i].sigma_rezi_PAS[it] = (double *) Malloc(sizeof(double) * NDIM, "IonsInitRead: sigma_rezi_PAS");
|
|---|
| 321 | }
|
|---|
| 322 | //I->I[i].IonFac = I->I[i].IonMass;
|
|---|
| [961b75] | 323 | I->I[i].IonMass *= Units2Electronmass; // in config given in amu not in "atomic units" (electron mass, m_e = 1)
|
|---|
| 324 | // proportional factor between electron and proton mass
|
|---|
| 325 | I->I[i].IonFac = Units2Electronmass/I->I[i].IonMass;
|
|---|
| [a0bcf1] | 326 | I->TotalMass += I->I[i].Max_IonsOfType*I->I[i].IonMass;
|
|---|
| [961b75] | 327 | sigma = sqrt(P->R.TargetTemp/I->I[i].IonMass);
|
|---|
| [a0bcf1] | 328 | I->I[i].ZFactor = -I->I[i].Z*4.*PI;
|
|---|
| 329 | I->I[i].alpha = (double *) Malloc(sizeof(double) * I->I[i].Max_IonsOfType, "IonsInitRead: alpha");
|
|---|
| 330 | I->I[i].R = (double *) Malloc(sizeof(double) * NDIM * I->I[i].Max_IonsOfType, "IonsInitRead: R");
|
|---|
| 331 | I->I[i].R_old = (double *) Malloc(sizeof(double) *NDIM* I->I[i].Max_IonsOfType, "IonsInitRead: R_old");
|
|---|
| 332 | I->I[i].R_old_old = (double *) Malloc(sizeof(double) *NDIM* I->I[i].Max_IonsOfType, "IonsInitRead: R_old_old");
|
|---|
| [961b75] | 333 | I->I[i].Constraints = (struct IonConstrained **) Malloc(sizeof(struct IonConstrained *)*I->I[i].Max_IonsOfType, "IonsInitRead: **Constraints");
|
|---|
| 334 | for (it=0;it<I->I[i].Max_IonsOfType;it++)
|
|---|
| 335 | I->I[i].Constraints[it] = NULL;
|
|---|
| [a0bcf1] | 336 | I->I[i].FIon = (double *) Malloc(sizeof(double) * NDIM * I->I[i].Max_IonsOfType, "IonsInitRead: FIon");
|
|---|
| 337 | I->I[i].FIon_old = (double *) Malloc(sizeof(double) * NDIM * I->I[i].Max_IonsOfType, "IonsInitRead: FIon_old");
|
|---|
| 338 | SetArrayToDouble0(I->I[i].FIon_old, NDIM*I->I[i].Max_IonsOfType);
|
|---|
| 339 | I->I[i].SearchDir = Malloc(sizeof(double) * NDIM * I->I[i].Max_IonsOfType, "IonsInitRead: SearchDir");
|
|---|
| 340 | I->I[i].GammaA = Malloc(sizeof(double) * I->I[i].Max_IonsOfType, "IonsInitRead: GammaA");
|
|---|
| 341 | SetArrayToDouble0(I->I[i].GammaA, I->I[i].Max_IonsOfType);
|
|---|
| 342 | I->I[i].U = (double *) Malloc(sizeof(double) *NDIM* I->I[i].Max_IonsOfType, "IonsInitRead: U");
|
|---|
| 343 | SetArrayToDouble0(I->I[i].U, NDIM*I->I[i].Max_IonsOfType);
|
|---|
| 344 | I->I[i].SFactor = NULL;
|
|---|
| 345 | I->I[i].IMT = (enum IonMoveType *) Malloc(sizeof(enum IonMoveType) * I->I[i].Max_IonsOfType, "IonsInitRead: IMT");
|
|---|
| 346 | I->I[i].FIonL = (double *) Malloc(sizeof(double) * NDIM * I->I[i].Max_IonsOfType, "IonsInitRead: FIonL");
|
|---|
| 347 | I->I[i].FIonNL = (double *) Malloc(sizeof(double) * NDIM * I->I[i].Max_IonsOfType, "IonsInitRead: FIonNL");
|
|---|
| [961b75] | 348 | I->I[i].FMagnetic = (double *) Malloc(sizeof(double) * NDIM * I->I[i].Max_IonsOfType, "IonsInitRead: FMagnetic");
|
|---|
| 349 | I->I[i].FConstraint = (double *) Malloc(sizeof(double) * NDIM * I->I[i].Max_IonsOfType, "IonsInitRead: FConstraint");
|
|---|
| [a0bcf1] | 350 | I->I[i].FEwald = (double *) Malloc(sizeof(double) * NDIM * I->I[i].Max_IonsOfType, "IonsInitRead: FEwald");
|
|---|
| 351 | I->I[i].alpha = (double *) Malloc(sizeof(double) * I->I[i].Max_IonsOfType, "IonsInitRead: alpha");
|
|---|
| [961b75] | 352 | // now parse ion coordination
|
|---|
| [a0bcf1] | 353 | for (j=0; j < I->I[i].Max_IonsOfType; j++) {
|
|---|
| 354 | I->I[i].alpha[j] = 2.;
|
|---|
| 355 | sprintf(name,"Ion_Type%i_%i",i+1,j+1);
|
|---|
| [961b75] | 356 | for (d=0; d <NDIM; d++) // transfor to old and old_old coordinates as well
|
|---|
| 357 | I->I[i].FMagnetic[d+j*NDIM] = 0.;
|
|---|
| 358 | // first ones are starting positions
|
|---|
| 359 | if ((count = ParseIonsCoordinatesAndVelocities(P, source, name, 1, relative, sigma, &I->I[i].R[NDIM*j], &I->I[i].U[NDIM*j], &IMT))) {
|
|---|
| 360 | for (d=0; d <NDIM; d++) // transfor to old and old_old coordinates as well
|
|---|
| 361 | I->I[i].R_old_old[d+NDIM*j] = I->I[i].R_old[d+NDIM*j] = I->I[i].R[d+NDIM*j];
|
|---|
| 362 | I->I[i].IMT[j] = IMT;
|
|---|
| 363 | if (P->Call.out[ReadOut]) fprintf(stderr, "(%i) R[%d,%d] = (%lg, %lg, %lg)\n", P->Par.me, i,j,I->I[i].R[0+NDIM*j],I->I[i].R[1+NDIM*j],I->I[i].R[2+NDIM*j]);
|
|---|
| 364 | if (P->Call.out[ReadOut]) fprintf(stderr, "(%i) U[%d,%d] = (%lg, %lg, %lg)\n", P->Par.me, i,j,I->I[i].U[0+NDIM*j],I->I[i].U[1+NDIM*j],I->I[i].U[2+NDIM*j]);
|
|---|
| 365 | if (P->Call.out[ReadOut]) fprintf(stderr, "(%i) IMT[%d,%d] = %i\n", P->Par.me, i,j,I->I[i].IMT[j]);
|
|---|
| [a0bcf1] | 366 | } else {
|
|---|
| [961b75] | 367 | Error(SomeError, "Could not find or parse coordinates of an ion");
|
|---|
| [a0bcf1] | 368 | }
|
|---|
| [961b75] | 369 | // following ones are constrained motions
|
|---|
| 370 | while (ParseIonsCoordinatesAndVelocities(P, source, name, ++count, relative, sigma, R, U, &IMT)) {
|
|---|
| 371 | if (P->Call.out[ReadOut]) fprintf(stderr, "(%i) Constrained read %d ...\n", P->Par.me, count-1);
|
|---|
| 372 | ptr = AddConstraintItem(&I->I[i], j); // allocate memory
|
|---|
| 373 | for (d=0; d< NDIM; d++) { // fill the constraint item
|
|---|
| 374 | ptr->R[d] = R[d];
|
|---|
| 375 | ptr->U[d] = U[d];
|
|---|
| [a0bcf1] | 376 | }
|
|---|
| [961b75] | 377 | ptr->IMT = IMT;
|
|---|
| 378 | ptr->step = count-1;
|
|---|
| 379 | if (P->Call.out[ReadOut]) fprintf(stderr, "(%i) R[%d,%d] = (%lg, %lg, %lg)\n", P->Par.me, i,j,ptr->R[0],ptr->R[1],ptr->R[2]);
|
|---|
| 380 | if (P->Call.out[ReadOut]) fprintf(stderr, "(%i) U[%d,%d] = (%lg, %lg, %lg)\n", P->Par.me, i,j,ptr->U[0],ptr->U[1],ptr->U[2]);
|
|---|
| 381 | if (P->Call.out[ReadOut]) fprintf(stderr, "(%i) IMT[%d,%d] = %i\n", P->Par.me, i,j,ptr->IMT);
|
|---|
| [a0bcf1] | 382 | }
|
|---|
| [961b75] | 383 | //if (P->Call.out[ReadOut])
|
|---|
| 384 | fprintf(stderr,"(%i) A total of %d additional constrained moves for ion (%d,%d) was parsed.\n", P->Par.me, count-2, i,j);
|
|---|
| [a0bcf1] | 385 |
|
|---|
| 386 | I->Max_TotalIons++;
|
|---|
| 387 | }
|
|---|
| 388 | }
|
|---|
| [32de28] | 389 | Free(free_name, "IonsInitRead: free_name");
|
|---|
| [a0bcf1] | 390 | I->Max_Max_IonsOfType = 0;
|
|---|
| 391 | for (i=0; i < I->Max_Types; i++)
|
|---|
| 392 | I->Max_Max_IonsOfType = MAX(I->Max_Max_IonsOfType, I->I[i].Max_IonsOfType);
|
|---|
| 393 | I->FTemp = (double *) Malloc(sizeof(double) * NDIM * I->Max_Max_IonsOfType, "IonsInitRead:");
|
|---|
| 394 | }
|
|---|
| 395 |
|
|---|
| [32de28] | 396 | /** Allocates memory for a IonConstrained item and adds to end of list.
|
|---|
| 397 | * \param *I IonType structure
|
|---|
| 398 | * \param ion which ion of the type
|
|---|
| 399 | * \return pointer to created item
|
|---|
| 400 | */
|
|---|
| 401 | struct IonConstrained * AddConstraintItem(struct IonType *I, int ion)
|
|---|
| 402 | {
|
|---|
| 403 | struct IonConstrained **ptr = &(I->Constraints[ion]);
|
|---|
| 404 | while (*ptr != NULL) { // advance to end of list
|
|---|
| 405 | ptr = &((*ptr)->next);
|
|---|
| 406 | }
|
|---|
| 407 | // allocated memory for structure and items within
|
|---|
| 408 | *ptr = (struct IonConstrained *) Malloc(sizeof(struct IonConstrained), "CreateConstraintItem: IonConstrained");
|
|---|
| 409 | (*ptr)->R = (double *) Malloc(sizeof(double)*NDIM, "AddConstraintItem: *R");
|
|---|
| 410 | (*ptr)->U = (double *) Malloc(sizeof(double)*NDIM, "AddConstraintItem: *U");
|
|---|
| 411 | (*ptr)->next = NULL;
|
|---|
| 412 | return *ptr;
|
|---|
| 413 | }
|
|---|
| 414 |
|
|---|
| 415 | /** Removes first of item of IonConstrained list for given \a ion.
|
|---|
| 416 | * Frees memory of items within and then
|
|---|
| 417 | * \param *I IonType structure
|
|---|
| 418 | * \param ion which ion of the type
|
|---|
| 419 | * \return 1 - success, 0 - failure (list is already empty, start pointer was NULL)
|
|---|
| 420 | */
|
|---|
| 421 | int RemoveConstraintItem(struct IonType *I, int ion)
|
|---|
| 422 | {
|
|---|
| 423 | struct IonConstrained **ptr = &(I->Constraints[ion]);
|
|---|
| 424 | struct IonConstrained *next_ptr;
|
|---|
| 425 | if (*ptr != NULL) {
|
|---|
| 426 | next_ptr = (*ptr)->next;
|
|---|
| 427 | Free((*ptr)->R, "RemoveConstraintItem: R");
|
|---|
| 428 | Free((*ptr)->U, "RemoveConstraintItem: U");
|
|---|
| 429 | Free((*ptr), "RemoveConstraintItem: IonConstrained structure");
|
|---|
| 430 | *ptr = next_ptr; // set start of list to next item (automatically is null if ptr was last item)
|
|---|
| 431 | return 1;
|
|---|
| 432 | } else {
|
|---|
| 433 | return 0;
|
|---|
| 434 | }
|
|---|
| 435 | }
|
|---|
| 436 |
|
|---|
| [a0bcf1] | 437 | /** Calculated Ewald force.
|
|---|
| [c5bdb23] | 438 | * The basic idea of the ewald method is to add a countercharge - here of gaussian form -
|
|---|
| 439 | * which shields the long-range charges making them effectively short-ranged, while the
|
|---|
| 440 | * countercharges themselves can be analytically (here by (hidden) Fouriertransform: it's
|
|---|
| 441 | * added to the electron density) integrated and withdrawn.
|
|---|
| [a0bcf1] | 442 | * \f[
|
|---|
| 443 | * E_{K-K} = \frac{1}{2} \sum_L \sum_{K=\{(I_s,I_a,J_s,J_a)|R_{I_s,I_a} - R_{J_s,J_a} - L\neq 0\}} \frac{Z_{I_s} Z_{J_s}} {|R_{I_s,I_a} - R_{J_s,J_a} - L|}
|
|---|
| 444 | * \textnormal{erfc} \Bigl ( \frac{|R_{I_s,I_a} - R_{J_s,J_a} - L|} {\sqrt{r_{I_s}^{Gauss}+r_{J_s}^{Gauss}}}\Bigr )
|
|---|
| 445 | * \qquad (4.10)
|
|---|
| 446 | * \f]
|
|---|
| 447 | * Calculates the core-to-core force between the ions of all super cells.
|
|---|
| 448 | * In order for this series to converge there must be a certain summation
|
|---|
| 449 | * applied, which is the ewald summation and which is nothing else but to move
|
|---|
| 450 | * in a circular motion from the current cell to the outside up to Ions::R_cut.
|
|---|
| 451 | * \param *P Problem at hand
|
|---|
| [c5bdb23] | 452 | * \param first if != 0 table mirror cells to take into (local!) account of ewald summation
|
|---|
| [a0bcf1] | 453 | */
|
|---|
| 454 | void CalculateEwald(struct Problem *P, int first)
|
|---|
| 455 | {
|
|---|
| [c5bdb23] | 456 | int r,is1,is2,ia1,ia2,ir,ir2,i,j,k,ActualVec,MaxVec,MaxCell,Is_Neighb_Cell,LocalActualVec; //,is
|
|---|
| [a0bcf1] | 457 | int MaxPar=P->Par.procs, ParMe=P->Par.me, MaxLocalVec, StartVec;
|
|---|
| 458 | double rcut2,r2,erre2,addesr,addpre,arg,fac,gkl,esr,fxx,repand;
|
|---|
| 459 | double R1[3];
|
|---|
| 460 | double R2[3];
|
|---|
| 461 | double R3[3];
|
|---|
| 462 | double t[3];
|
|---|
| 463 | struct Lattice *L = &P->Lat;
|
|---|
| 464 | struct PseudoPot *PP = &P->PP;
|
|---|
| 465 | struct Ions *I = &P->Ion;
|
|---|
| [c5bdb23] | 466 | if (first) {
|
|---|
| 467 | SpeedMeasure(P, InitSimTime, StopTimeDo);
|
|---|
| 468 | rcut2 = I->R_cut*I->R_cut;
|
|---|
| 469 | ActualVec = 0; // counts number of cells taken into account
|
|---|
| 470 | MaxCell = (int)(5.*I->R_cut/pow(L->Volume,1./3.)); // the number of cells in each direction is prop. to axis length over cutoff
|
|---|
| 471 | if (MaxCell < 2) MaxCell = 2; // take at least 2
|
|---|
| 472 | for (i = -MaxCell; i <= MaxCell; i++) {
|
|---|
| 473 | for (j = -MaxCell; j <= MaxCell; j++) {
|
|---|
| 474 | for (k = -MaxCell; k <= MaxCell; k++) {
|
|---|
| 475 | r2 = 0.;
|
|---|
| 476 | for (ir=0; ir <NDIM; ir++) { // t is the offset vector pointing to distant cell (only diagonal entries)
|
|---|
| 477 | t[ir] = i*L->RealBasis[0*NDIM+ir]+j*L->RealBasis[1*NDIM+ir]+k*L->RealBasis[2*NDIM+ir];
|
|---|
| 478 | r2 = t[ir]*t[ir]; // is squared distance to other cell
|
|---|
| 479 | }
|
|---|
| 480 | Is_Neighb_Cell = ((abs(i)<=1) && (abs(j)<=1) && (abs(k)<=1)); // whether it's directly adjacent
|
|---|
| 481 | if ((r2 <= rcut2) || Is_Neighb_Cell) { // if it's either adjacent or within cutoff, count it in
|
|---|
| 482 | ActualVec++;
|
|---|
| 483 | }
|
|---|
| 484 | }
|
|---|
| [a0bcf1] | 485 | }
|
|---|
| [c5bdb23] | 486 | }
|
|---|
| 487 | MaxVec = ActualVec; // store number of cells
|
|---|
| 488 | I->MaxVec = ActualVec;
|
|---|
| 489 | MaxLocalVec = MaxVec / MaxPar; // share among processes
|
|---|
| 490 | StartVec = ParMe*MaxLocalVec; // for this process start at its patch
|
|---|
| 491 | r = MaxVec % MaxPar; // rest not fitting...
|
|---|
| 492 | if (ParMe >= r) { // ones who don't get extra cells have to shift their start vector
|
|---|
| 493 | StartVec += r;
|
|---|
| 494 | } else { // others get extra cell and have to shift also by a bit
|
|---|
| 495 | StartVec += ParMe;
|
|---|
| 496 | MaxLocalVec++;
|
|---|
| 497 | }
|
|---|
| 498 | I->MaxLocalVec = MaxLocalVec;
|
|---|
| 499 | LocalActualVec = ActualVec = 0;
|
|---|
| 500 | // now go through the same loop again and note down every offset vector for cells in this process' patch
|
|---|
| 501 | if (MaxLocalVec != 0) {
|
|---|
| [a0bcf1] | 502 | I->RLatticeVec = (double *)
|
|---|
| [c5bdb23] | 503 | Realloc(I->RLatticeVec, sizeof(double)*NDIM*MaxLocalVec, "CalculateEwald: I->RLatticeVec");
|
|---|
| 504 | } else {
|
|---|
| 505 | fprintf(stderr,"(%i) Warning in Ewald summation: Got MaxLocalVec == 0\n", P->Par.me);
|
|---|
| 506 | }
|
|---|
| 507 | for (i = -MaxCell; i <= MaxCell; i++) {
|
|---|
| 508 | for (j = -MaxCell; j <= MaxCell; j++) {
|
|---|
| 509 | for (k = -MaxCell; k <= MaxCell; k++) {
|
|---|
| 510 | r2 = 0.;
|
|---|
| 511 | for (ir=0; ir <NDIM; ir++) { // create offset vector from diagonal entries
|
|---|
| 512 | t[ir] = i*L->RealBasis[0*NDIM+ir]+j*L->RealBasis[1*NDIM+ir]+k*L->RealBasis[2*NDIM+ir];
|
|---|
| 513 | r2 = t[ir]*t[ir];
|
|---|
| 514 | }
|
|---|
| 515 | Is_Neighb_Cell = ((abs(i)<=1) && (abs(j)<=1) && (abs(k)<=1));
|
|---|
| 516 | if ((r2 <= rcut2) || Is_Neighb_Cell) { // if adjacent or within cutoff ...
|
|---|
| 517 | if (ActualVec >= StartVec && ActualVec < StartVec+MaxLocalVec) { // ... and belongs to patch, store 'em
|
|---|
| 518 | I->RLatticeVec[0+NDIM*LocalActualVec] = t[0];
|
|---|
| 519 | I->RLatticeVec[1+NDIM*LocalActualVec] = t[1];
|
|---|
| 520 | I->RLatticeVec[2+NDIM*LocalActualVec] = t[2];
|
|---|
| 521 | LocalActualVec++;
|
|---|
| 522 | }
|
|---|
| 523 | ActualVec++;
|
|---|
| 524 | }
|
|---|
| 525 | }
|
|---|
| [a0bcf1] | 526 | }
|
|---|
| 527 | }
|
|---|
| [c5bdb23] | 528 | SpeedMeasure(P, InitSimTime, StartTimeDo);
|
|---|
| 529 | }
|
|---|
| 530 | SpeedMeasure(P, EwaldTime, StartTimeDo);
|
|---|
| 531 |
|
|---|
| 532 | // first set everything to zero
|
|---|
| 533 | esr = 0.0;
|
|---|
| 534 | //for (is1=0;is1 < I->Max_Types; is1++)
|
|---|
| 535 | // then take each ion of each type once
|
|---|
| 536 | for (is1=0;is1 < I->Max_Types; is1++) {
|
|---|
| 537 | // reset FTemp for IonType
|
|---|
| 538 | for (ia1=0;ia1 < I->I[is1].Max_IonsOfType; ia1++)
|
|---|
| 539 | for (i=0; i< NDIM; i++)
|
|---|
| 540 | I->FTemp[i+NDIM*(ia1)] = 0.;
|
|---|
| 541 | // then sum for each on of the IonType
|
|---|
| 542 | for (ia1=0;ia1 < I->I[is1].Max_IonsOfType; ia1++) {
|
|---|
| 543 | for (i=0;i<NDIM;i++) {
|
|---|
| 544 | R1[i]=I->I[is1].R[i+NDIM*ia1]; // R1 is local coordinate of current first ion
|
|---|
| 545 | }
|
|---|
| 546 | for (ir2=0;ir2<I->MaxLocalVec;ir2++) { // for every cell of local patch (each with the same atoms)
|
|---|
| 547 | for (is2=0;is2 < I->Max_Types; is2++) { // and for each ion of each type a second time
|
|---|
| 548 | // gkl = 1./(sqrt(r_is1,gauss^2 + r_is2,gauss^2))
|
|---|
| 549 | gkl=1./sqrt(I->I[is1].rgauss*I->I[is1].rgauss + I->I[is2].rgauss*I->I[is2].rgauss);
|
|---|
| 550 | for (ia2=0;ia2<I->I[is2].Max_IonsOfType; ia2++) {
|
|---|
| 551 | for (i=0;i<3;i++) {
|
|---|
| 552 | R2[i]=I->RLatticeVec[i+NDIM*ir2]+I->I[is2].R[i+NDIM*ia2]; // R2 is coordinate of current second ion in the distant cell!
|
|---|
| 553 | }
|
|---|
| 554 |
|
|---|
| 555 | erre2=0.0;
|
|---|
| 556 | for (i=0;i<NDIM;i++) {
|
|---|
| 557 | R3[i] = R1[i]-R2[i]; // R3 = R_is1,ia1 - R_is2,ia2 - L
|
|---|
| 558 | erre2+=R3[i]*R3[i]; // erre2 = |R_is1,ia1 - R_is2,ia2 - L|^2
|
|---|
| 559 | }
|
|---|
| 560 | if (erre2 > MYEPSILON) { // appears as one over ... in fac, thus check if zero
|
|---|
| 561 | arg = sqrt(erre2); // arg = |R_is1,ia1 - R_is2,ia2 - L|
|
|---|
| 562 | fac=PP->zval[is1]*PP->zval[is2]/arg*0.5; // fac = -1/2 * Z_is1 * Z_is2 / arg
|
|---|
| 563 |
|
|---|
| 564 | arg *= gkl; // arg = |R_is1,ia1 - R_is2,ia2 - L|/(sqrt(r_is1,gauss^2 + r_is2,gauss^2))
|
|---|
| 565 | addesr = derf(arg);
|
|---|
| 566 | addesr = (1.-addesr)*fac; // complementary error function: addesr = 1/2*erfc(arg)/|R_is1,ia1 - R_is2,ia2 - L|
|
|---|
| 567 | esr += addesr;
|
|---|
| 568 | addpre=exp(-arg*arg)*gkl;
|
|---|
| 569 | addpre = PP->fac1sqrtPI*PP->zval[is1]*PP->zval[is2]*addpre; // addpre = exp(-|R_is1,ia1 - R_is2,ia2 - L|^2/(r_is1,gauss^2 + r_is2,gauss^2))/(sqrt(pi)*sqrt(r_is1,gauss^2 + r_is2,gauss^2))
|
|---|
| 570 | repand = (addesr+addpre)/erre2;
|
|---|
| 571 | // if ((fabs(repand) > MYEPSILON)) {
|
|---|
| 572 | // fprintf(stderr,"(%i) Ewald[is%d,ia%d/is%d,ia%d,(%d)]: Z1 %lg, Z2 %lg\tpre %lg\t esr %lg\t fac %lg\t arg %lg\tR3 (%lg,%lg,%lg)\n", P->Par.me, is1, ia1, is2, ia2, ir2, PP->zval[is1], PP->zval[is2],addpre,addesr, fac, arg, R3[0], R3[1], R3[2]);
|
|---|
| 573 | // }
|
|---|
| 574 | for (i=0;i<NDIM;i++) {
|
|---|
| 575 | fxx=repand*R3[i];
|
|---|
| 576 | // if ((fabs(repand) > MYEPSILON)) {
|
|---|
| 577 | // fprintf(stderr,"%i %i %i %i %i %g\n",is1+1,ia1+1,is2+1,ia2+1,i+1,fxx);
|
|---|
| 578 | // }
|
|---|
| 579 | I->FTemp[i+NDIM*(ia1)] += fxx*2.0; // actio = reactio?
|
|---|
| 580 | //I->FTemp[i+NDIM*(ia2)] -= fxx;
|
|---|
| 581 | }
|
|---|
| 582 | }
|
|---|
| 583 | }
|
|---|
| 584 | }
|
|---|
| [a0bcf1] | 585 | }
|
|---|
| 586 | }
|
|---|
| [c5bdb23] | 587 | MPI_Allreduce (I->FTemp , I->I[is1].FEwald, NDIM*I->I[is1].Max_IonsOfType, MPI_DOUBLE, MPI_SUM, P->Par.comm);
|
|---|
| [a0bcf1] | 588 | }
|
|---|
| 589 | SpeedMeasure(P, EwaldTime, StopTimeDo);
|
|---|
| [c5bdb23] | 590 | //for (is=0;is < I->Max_Types; is++) {
|
|---|
| 591 | //}
|
|---|
| [a0bcf1] | 592 | MPI_Allreduce ( &esr, &L->E->AllTotalIonsEnergy[EwaldEnergy], 1, MPI_DOUBLE, MPI_SUM, P->Par.comm);
|
|---|
| [c5bdb23] | 593 | // fprintf(stderr,"\n");
|
|---|
| [a0bcf1] | 594 | }
|
|---|
| 595 |
|
|---|
| 596 | /** Sum up all forces acting on Ions.
|
|---|
| [69eca8] | 597 | * IonType::FIon = IonType::FEwald + IonType::FIonL + IonType::FIonNL + IonType::FMagnetic for all
|
|---|
| [a0bcf1] | 598 | * dimensions #NDIM and each Ion of IonType
|
|---|
| 599 | * \param *P Problem at hand
|
|---|
| 600 | */
|
|---|
| 601 | void CalculateIonForce(struct Problem *P)
|
|---|
| 602 | {
|
|---|
| 603 | struct Ions *I = &P->Ion;
|
|---|
| 604 | int i,j,d;
|
|---|
| 605 | for (i=0; i < I->Max_Types; i++)
|
|---|
| 606 | for (j=0; j < I->I[i].Max_IonsOfType; j++)
|
|---|
| 607 | for (d=0; d<NDIM;d++)
|
|---|
| [69eca8] | 608 | I->I[i].FIon[d+j*NDIM] = (I->I[i].FEwald[d+j*NDIM] + I->I[i].FIonL[d+j*NDIM] + I->I[i].FIonNL[d+j*NDIM] + I->I[i].FMagnetic[d+j*NDIM]);
|
|---|
| [a0bcf1] | 609 | }
|
|---|
| 610 |
|
|---|
| 611 | /** Write Forces to FileData::ForcesFile.
|
|---|
| 612 | * goes through all Ions per IonType and each dimension of #NDIM and prints in one line:
|
|---|
| 613 | * Position IonType::R, total force Ion IonType::FIon, local force IonType::FIonL,
|
|---|
| 614 | * non-local IonType::FIonNL, ewald force IonType::FEwald
|
|---|
| 615 | * \param *P Problem at hand
|
|---|
| 616 | */
|
|---|
| 617 | void OutputIonForce(struct Problem *P)
|
|---|
| 618 | {
|
|---|
| 619 | struct RunStruct *R = &P->R;
|
|---|
| 620 | struct FileData *F = &P->Files;
|
|---|
| 621 | struct Ions *I = &P->Ion;
|
|---|
| 622 | FILE *fout = NULL;
|
|---|
| 623 | int i,j;
|
|---|
| 624 | if (F->MeOutMes != 1) return;
|
|---|
| 625 | fout = F->ForcesFile;
|
|---|
| 626 | for (i=0; i < I->Max_Types; i++)
|
|---|
| 627 | for (j=0; j < I->I[i].Max_IonsOfType; j++)
|
|---|
| [8dd187] | 628 | fprintf(fout, "%i\t%i\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\t%e\n", i,j,
|
|---|
| [a0bcf1] | 629 | I->I[i].R[0+j*NDIM], I->I[i].R[1+j*NDIM], I->I[i].R[2+j*NDIM],
|
|---|
| 630 | I->I[i].FIon[0+j*NDIM], I->I[i].FIon[1+j*NDIM], I->I[i].FIon[2+j*NDIM],
|
|---|
| 631 | I->I[i].FIonL[0+j*NDIM], I->I[i].FIonL[1+j*NDIM], I->I[i].FIonL[2+j*NDIM],
|
|---|
| 632 | I->I[i].FIonNL[0+j*NDIM], I->I[i].FIonNL[1+j*NDIM], I->I[i].FIonNL[2+j*NDIM],
|
|---|
| [8dd187] | 633 | I->I[i].FMagnetic[0+j*NDIM], I->I[i].FMagnetic[1+j*NDIM], I->I[i].FMagnetic[2+j*NDIM],
|
|---|
| [a0bcf1] | 634 | I->I[i].FEwald[0+j*NDIM], I->I[i].FEwald[1+j*NDIM], I->I[i].FEwald[2+j*NDIM]);
|
|---|
| 635 | fflush(fout);
|
|---|
| 636 | fprintf(fout, "MeanForce:\t%e\n", R->MeanForce[0]);
|
|---|
| 637 | }
|
|---|
| 638 |
|
|---|
| [b74e5e] | 639 | /** Reads Forces from CallOptions::ForcesFile.
|
|---|
| 640 | * goes through all Ions per IonType and each dimension of #NDIM and parses in one line:
|
|---|
| 641 | * Position IonType::R, total force Ion IonType::FIon, local force IonType::FIonL,
|
|---|
| 642 | * non-local IonType::FIonNL, ewald force IonType::FEwald
|
|---|
| 643 | * \param *P Problem at hand
|
|---|
| 644 | */
|
|---|
| 645 | void ParseIonForce(struct Problem *P)
|
|---|
| 646 | {
|
|---|
| 647 | //struct RunStruct *Run = &P->R;
|
|---|
| 648 | struct Ions *I = &P->Ion;
|
|---|
| 649 | FILE *finput = NULL;
|
|---|
| 650 | char line[1024];
|
|---|
| 651 | int i,j;
|
|---|
| 652 | double i1, j1;
|
|---|
| 653 | double R[NDIM];
|
|---|
| 654 |
|
|---|
| 655 | if (P->Call.ForcesFile == NULL) return;
|
|---|
| 656 | finput = fopen(P->Call.ForcesFile, "r");
|
|---|
| 657 | //fprintf(stderr, "File pointer says %p\n", finput);
|
|---|
| 658 | if (finput == NULL)
|
|---|
| 659 | Error(InitReading, "ForcesFile does not exist.");
|
|---|
| 660 | fscanf(finput, "%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\t%s\n", line,line,line,line,line,line,line,line,line,line,line,line,line,line,line,line,line,line,line,line); // skip header line
|
|---|
| 661 | for (i=0; i < I->Max_Types; i++)
|
|---|
| 662 | for (j=0; j < I->I[i].Max_IonsOfType; j++)
|
|---|
| 663 | if (feof(finput)) {
|
|---|
| 664 | Error(InitReading, "ForcesFile does not contain enough lines.");
|
|---|
| 665 | } else {
|
|---|
| 666 | //fscanf(finput, "%s\n", line);
|
|---|
| 667 | //if (strlen(line) < 100) {
|
|---|
| 668 | //fprintf(stderr, "i %d, j %d, length of line %ld, line %s.\n", i,j, strlen(line), line);
|
|---|
| 669 | //Error(InitReading, "ForcesFile does not contain enough lines or line is faulty.");
|
|---|
| 670 | //} else
|
|---|
| 671 | //fprintf(stderr, "Parsed line: '%s'\n", line);
|
|---|
| 672 | fscanf(finput, "%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\t%le\n", &i1,&j1,
|
|---|
| 673 | &R[0], &R[1], &R[2],
|
|---|
| 674 | &I->I[i].FIon[0+j*NDIM], &I->I[i].FIon[1+j*NDIM], &I->I[i].FIon[2+j*NDIM],
|
|---|
| 675 | &I->I[i].FIonL[0+j*NDIM], &I->I[i].FIonL[1+j*NDIM], &I->I[i].FIonL[2+j*NDIM],
|
|---|
| 676 | &I->I[i].FIonNL[0+j*NDIM], &I->I[i].FIonNL[1+j*NDIM], &I->I[i].FIonNL[2+j*NDIM],
|
|---|
| 677 | &I->I[i].FMagnetic[0+j*NDIM], &I->I[i].FMagnetic[1+j*NDIM], &I->I[i].FMagnetic[2+j*NDIM],
|
|---|
| 678 | &I->I[i].FEwald[0+j*NDIM], &I->I[i].FEwald[1+j*NDIM], &I->I[i].FEwald[2+j*NDIM]);
|
|---|
| 679 | if ((i != (int)i1) || (j != (int)j1))
|
|---|
| 680 | Error(InitReading, "Line does not match to desired ion.");
|
|---|
| 681 | }
|
|---|
| 682 | fclose(finput);
|
|---|
| 683 | //fprintf(stderr, "MeanForce:\t%e\n", Run->MeanForce[0]);
|
|---|
| 684 | };
|
|---|
| 685 |
|
|---|
| [a0bcf1] | 686 | /** Frees memory IonType.
|
|---|
| 687 | * All pointers from IonType and the pointer on it are free'd
|
|---|
| 688 | * \param *I Ions structure to be free'd
|
|---|
| 689 | * \sa IonsInitRead where memory is allocated
|
|---|
| 690 | */
|
|---|
| 691 | void RemoveIonsRead(struct Ions *I)
|
|---|
| 692 | {
|
|---|
| 693 | int i, it;
|
|---|
| 694 | for (i=0; i < I->Max_Types; i++) {
|
|---|
| [32de28] | 695 | Free(I->I[i].Name, "RemoveIonsRead: I->I[i].Name");
|
|---|
| 696 | Free(I->I[i].Symbol, "RemoveIonsRead: I->I[i].Symbol");
|
|---|
| [a0bcf1] | 697 | for (it=0;it<I->I[i].Max_IonsOfType;it++) {
|
|---|
| [32de28] | 698 | Free(I->I[i].moment[it], "RemoveIonsRead: I->I[i].moment[it]");
|
|---|
| 699 | Free(I->I[i].sigma[it], "RemoveIonsRead: I->I[i].sigma[it]");
|
|---|
| 700 | Free(I->I[i].sigma_rezi[it], "RemoveIonsRead: I->I[i].sigma_rezi[it]");
|
|---|
| 701 | Free(I->I[i].sigma_PAS[it], "RemoveIonsRead: I->I[i].sigma_PAS[it]");
|
|---|
| 702 | Free(I->I[i].sigma_rezi_PAS[it], "RemoveIonsRead: I->I[i].sigma_rezi_PAS[it]");
|
|---|
| 703 | while (RemoveConstraintItem(&I->I[i], it)); // empty the constrained item list
|
|---|
| [a0bcf1] | 704 | }
|
|---|
| [32de28] | 705 | Free(I->I[i].sigma, "RemoveIonsRead: I->I[i].sigma");
|
|---|
| 706 | Free(I->I[i].sigma_rezi, "RemoveIonsRead: I->I[i].sigma_rezi");
|
|---|
| 707 | Free(I->I[i].sigma_PAS, "RemoveIonsRead: I->I[i].sigma_PAS");
|
|---|
| 708 | Free(I->I[i].sigma_rezi_PAS, "RemoveIonsRead: I->I[i].sigma_rezi_PAS");
|
|---|
| 709 | Free(I->I[i].R, "RemoveIonsRead: I->I[i].R");
|
|---|
| 710 | Free(I->I[i].R_old, "RemoveIonsRead: I->I[i].R_old");
|
|---|
| 711 | Free(I->I[i].R_old_old, "RemoveIonsRead: I->I[i].R_old_old");
|
|---|
| 712 | Free(I->I[i].Constraints, "RemoveIonsRead: I->I[i].Constraints");
|
|---|
| 713 | Free(I->I[i].FIon, "RemoveIonsRead: I->I[i].FIon");
|
|---|
| 714 | Free(I->I[i].FIon_old, "RemoveIonsRead: I->I[i].FIon_old");
|
|---|
| 715 | Free(I->I[i].SearchDir, "RemoveIonsRead: I->I[i].SearchDir");
|
|---|
| 716 | Free(I->I[i].GammaA, "RemoveIonsRead: I->I[i].GammaA");
|
|---|
| 717 | Free(I->I[i].FIonL, "RemoveIonsRead: I->I[i].FIonL");
|
|---|
| 718 | Free(I->I[i].FIonNL, "RemoveIonsRead: I->I[i].FIonNL");
|
|---|
| 719 | Free(I->I[i].FMagnetic, "RemoveIonsRead: I->I[i].FMagnetic");
|
|---|
| 720 | Free(I->I[i].U, "RemoveIonsRead: I->I[i].U");
|
|---|
| 721 | Free(I->I[i].SFactor, "RemoveIonsRead: I->I[i].SFactor");
|
|---|
| 722 | Free(I->I[i].IMT, "RemoveIonsRead: I->I[i].IMT");
|
|---|
| 723 | Free(I->I[i].FEwald, "RemoveIonsRead: I->I[i].FEwald");
|
|---|
| 724 | Free(I->I[i].alpha, "RemoveIonsRead: I->I[i].alpha");
|
|---|
| [a0bcf1] | 725 | }
|
|---|
| 726 | if (I->R_cut == 0.0)
|
|---|
| [32de28] | 727 | Free(I->RLatticeVec, "RemoveIonsRead: I->RLatticeVec");
|
|---|
| 728 | Free(I->FTemp, "RemoveIonsRead: I->FTemp");
|
|---|
| 729 | Free(I->I, "RemoveIonsRead: I->I");
|
|---|
| [a0bcf1] | 730 | }
|
|---|
| 731 |
|
|---|
| 732 | /** Shifts center of gravity of ion forces IonType::FIon.
|
|---|
| [9a2c8c] | 733 | * First sums up all forces of movable ions for net center of gravity force,
|
|---|
| 734 | * then reduces each force by this temp over Ions#Max_TotalIons, so that all in all
|
|---|
| [a0bcf1] | 735 | * the net force is 0.
|
|---|
| 736 | * \param *P Problem at hand
|
|---|
| [9a2c8c] | 737 | * \sa CorrectVelocity()
|
|---|
| [a0bcf1] | 738 | */
|
|---|
| 739 | void CorrectForces(struct Problem *P)
|
|---|
| 740 | {
|
|---|
| 741 | struct Ions *I = &P->Ion;
|
|---|
| 742 | double *FIon;
|
|---|
| 743 | double FTemp[NDIM];
|
|---|
| [9a2c8c] | 744 | int is, ia, d, NumIons = 0;
|
|---|
| [a0bcf1] | 745 | for (d=0; d<NDIM; d++)
|
|---|
| [9a2c8c] | 746 | FTemp[d] = 0.;
|
|---|
| [a0bcf1] | 747 | for (is=0; is < I->Max_Types; is++) { // calculate force temperature for each type ...
|
|---|
| 748 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) { // .. and each ion of this type ...
|
|---|
| 749 | FIon = &I->I[is].FIon[NDIM*ia];
|
|---|
| [9a2c8c] | 750 | if (I->I[is].IMT[ia] == MoveIon) { // .. if it's influenced
|
|---|
| 751 | NumIons++;
|
|---|
| [a0bcf1] | 752 | for (d=0; d<NDIM; d++)
|
|---|
| 753 | FTemp[d] += FIon[d];
|
|---|
| 754 | }
|
|---|
| 755 | }
|
|---|
| 756 | }
|
|---|
| [9a2c8c] | 757 | //if (P->Files.MeOutMes != 1) return;
|
|---|
| 758 | // fprintf(stderr, "TotalForce before: %e\t %e\t %e\n", FTemp[0], FTemp[1], FTemp[2]);
|
|---|
| [a0bcf1] | 759 | for (is=0; is < I->Max_Types; is++) { // and then reduce each by this value
|
|---|
| 760 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
|---|
| 761 | FIon = &I->I[is].FIon[NDIM*ia];
|
|---|
| 762 | if (I->I[is].IMT[ia] == MoveIon) {
|
|---|
| 763 | for (d=0; d<NDIM; d++)
|
|---|
| [9a2c8c] | 764 | FIon[d] -= FTemp[d]/NumIons;
|
|---|
| [a0bcf1] | 765 | }
|
|---|
| 766 | }
|
|---|
| 767 | }
|
|---|
| [9a2c8c] | 768 | /* for (d=0; d<NDIM; d++)
|
|---|
| 769 | FTemp[d] = 0.;
|
|---|
| 770 | for (is=0; is < I->Max_Types; is++) { // calculate force temperature for each type ...
|
|---|
| 771 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) { // .. and each ion of this type ...
|
|---|
| 772 | FIon = &I->I[is].FIon[NDIM*ia];
|
|---|
| 773 | if (I->I[is].IMT[ia] == MoveIon) { // .. if it's influenced
|
|---|
| 774 | for (d=0; d<NDIM; d++)
|
|---|
| 775 | FTemp[d] += FIon[d];
|
|---|
| 776 | }
|
|---|
| 777 | }
|
|---|
| 778 | }
|
|---|
| 779 | if (P->Files.MeOutMes != 1) return;
|
|---|
| 780 | fprintf(stderr, "TotalForce after: %e\t %e\t %e\n", FTemp[0], FTemp[1], FTemp[2]);
|
|---|
| 781 | */
|
|---|
| [a0bcf1] | 782 | }
|
|---|
| 783 |
|
|---|
| 784 |
|
|---|
| 785 | /** Shifts center of gravity of ion velocities (rather momentums).
|
|---|
| 786 | * First sums up ion speed U times their IonMass (summed up for each
|
|---|
| [4782e78] | 787 | * dimension) in temp, then reduces the velocity by temp per Ions::TotalMass (so here
|
|---|
| [a0bcf1] | 788 | * total number of Ions is included)
|
|---|
| 789 | * \param *P Problem at hand
|
|---|
| [4782e78] | 790 | * \sa CorrectForces()
|
|---|
| [a0bcf1] | 791 | */
|
|---|
| 792 | void CorrectVelocity(struct Problem *P)
|
|---|
| 793 | {
|
|---|
| 794 | struct Ions *I = &P->Ion;
|
|---|
| 795 | double *U;
|
|---|
| 796 | double UTemp[NDIM];
|
|---|
| 797 | int is, ia, d;
|
|---|
| 798 | for (d=0; d<NDIM; d++)
|
|---|
| 799 | UTemp[d] = 0;
|
|---|
| 800 | for (is=0; is < I->Max_Types; is++) { // all types ...
|
|---|
| 801 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) { // ... all ions per type ...
|
|---|
| 802 | U = &I->I[is].U[NDIM*ia];
|
|---|
| [4782e78] | 803 | //if (I->I[is].IMT[ia] == MoveIon) { // ... even FixedIon moves, only not by other's forces
|
|---|
| [a0bcf1] | 804 | for (d=0; d<NDIM; d++)
|
|---|
| 805 | UTemp[d] += I->I[is].IonMass*U[d];
|
|---|
| [4782e78] | 806 | //}
|
|---|
| [a0bcf1] | 807 | }
|
|---|
| 808 | }
|
|---|
| 809 | for (is=0; is < I->Max_Types; is++) {
|
|---|
| 810 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
|---|
| 811 | U = &I->I[is].U[NDIM*ia];
|
|---|
| [4782e78] | 812 | //if (I->I[is].IMT[ia] == MoveIon) {
|
|---|
| [a0bcf1] | 813 | for (d=0; d<NDIM; d++)
|
|---|
| 814 | U[d] -= UTemp[d]/I->TotalMass; // shift by mean velocity over mass and number of ions
|
|---|
| [4782e78] | 815 | //}
|
|---|
| [a0bcf1] | 816 | }
|
|---|
| 817 | }
|
|---|
| 818 | }
|
|---|
| 819 |
|
|---|
| 820 | /** Moves ions according to calculated force.
|
|---|
| 821 | * Goes through each type of IonType and each ion therein, takes mass and
|
|---|
| 822 | * Newton to move each ion to new coordinates IonType::R, while remembering the last
|
|---|
| 823 | * two coordinates and the last force of the ion. Coordinates R are being
|
|---|
| 824 | * transformed to inverted base, shifted by +-1.0 and back-transformed to
|
|---|
| 825 | * real base.
|
|---|
| 826 | * \param *P Problem at hand
|
|---|
| 827 | * \sa CalculateIonForce
|
|---|
| 828 | * \sa UpdateIonsU
|
|---|
| 829 | * \warning U is not changed here, only used to move the ion.
|
|---|
| 830 | */
|
|---|
| 831 | void UpdateIonsR(struct Problem *P)
|
|---|
| 832 | {
|
|---|
| 833 | struct Ions *I = &P->Ion;
|
|---|
| 834 | int is, ia, d;
|
|---|
| 835 | double *R, *R_old, *R_old_old, *FIon, *FIon_old, *U, *FIonL, *FIonNL, *FEwald;
|
|---|
| 836 | double IonMass, a;
|
|---|
| 837 | double sR[NDIM], sRold[NDIM], sRoldold[NDIM];
|
|---|
| 838 | const int delta_t = P->R.delta_t;
|
|---|
| 839 | for (is=0; is < I->Max_Types; is++) { // go through all types ...
|
|---|
| 840 | IonMass = I->I[is].IonMass;
|
|---|
| 841 | if (IonMass < MYEPSILON) fprintf(stderr,"UpdateIonsR: IonMass = %lg",IonMass);
|
|---|
| 842 | a = delta_t*0.5/IonMass; // F/m = a and thus: s = 0.5 * F/m * t^2 + v * t =: t * (F * a + v)
|
|---|
| 843 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) { // .. and each ion of the type
|
|---|
| 844 | FIon = &I->I[is].FIon[NDIM*ia];
|
|---|
| 845 | FIonL = &I->I[is].FIonL[NDIM*ia];
|
|---|
| 846 | FIonNL = &I->I[is].FIonNL[NDIM*ia];
|
|---|
| 847 | FEwald = &I->I[is].FEwald[NDIM*ia];
|
|---|
| 848 | FIon_old = &I->I[is].FIon_old[NDIM*ia];
|
|---|
| 849 | U = &I->I[is].U[NDIM*ia];
|
|---|
| 850 | R = &I->I[is].R[NDIM*ia];
|
|---|
| 851 | R_old = &I->I[is].R_old[NDIM*ia];
|
|---|
| 852 | R_old_old = &I->I[is].R_old_old[NDIM*ia];
|
|---|
| [6590194] | 853 | //if (I->I[is].IMT[ia] == MoveIon) { // even FixedIon moves, only not by other's forces
|
|---|
| [a0bcf1] | 854 | for (d=0; d<NDIM; d++) {
|
|---|
| 855 | R_old_old[d] = R_old[d]; // shift old values
|
|---|
| 856 | R_old[d] = R[d];
|
|---|
| [6590194] | 857 | R[d] += delta_t*(U[d]+a*FIon[d]); // F = m * a and s = 0.5 * a * t^2
|
|---|
| [a0bcf1] | 858 | FIon_old[d] = FIon[d]; // store old force
|
|---|
| [6590194] | 859 | // FIon[d] = 0.; // zero all as a sign that's been moved
|
|---|
| 860 | // FIonL[d] = 0.;
|
|---|
| 861 | // FIonNL[d] = 0.;
|
|---|
| 862 | // FEwald[d] = 0.;
|
|---|
| 863 | // don't reset forces here anymore, as we OutputVis() after this and before next CalculateForce() (rendering them null in the cute graphs)
|
|---|
| [a0bcf1] | 864 | }
|
|---|
| [6590194] | 865 | if (I->I[is].Constraints[ia] != NULL) { // override current resulting R if constrained is given
|
|---|
| 866 | fprintf(stderr, "(%i) Using constrained coordinates R of step %d on ion (%i,%i)\n", P->Par.me, I->I[is].Constraints[ia]->step, is, ia);
|
|---|
| 867 | if ((I->I[is].Constraints[ia]->step != P->R.OuterStep) && (P->Par.me == 0))
|
|---|
| 868 | fprintf(stderr, "(%i) WARNING: constraint step %d does not match Outerstep %d for ion (%d,%d)\n", P->Par.me, I->I[is].Constraints[ia]->step, P->R.OuterStep, is, ia);
|
|---|
| 869 | for (d=0; d<NDIM; d++)
|
|---|
| 870 | R[d] = I->I[is].Constraints[ia]->R[d];
|
|---|
| 871 | I->I[is].IMT[ia] = I->I[is].Constraints[ia]->IMT;
|
|---|
| 872 | //RemoveConstraintItem(&I->I[is],ia); // and remove this first item from the list: // is done now at UpdateIonsU()
|
|---|
| 873 | }
|
|---|
| [a0bcf1] | 874 | RMat33Vec3(sR, P->Lat.InvBasis, R);
|
|---|
| 875 | RMat33Vec3(sRold, P->Lat.InvBasis, R_old);
|
|---|
| 876 | RMat33Vec3(sRoldold, P->Lat.InvBasis, R_old_old);
|
|---|
| 877 | for (d=0; d <NDIM; d++) {
|
|---|
| 878 | while (sR[d] < 0) {
|
|---|
| 879 | sR[d] += 1.0;
|
|---|
| 880 | sRold[d] += 1.0;
|
|---|
| 881 | sRoldold[d] += 1.0;
|
|---|
| 882 | }
|
|---|
| 883 | while (sR[d] >= 1.0) {
|
|---|
| 884 | sR[d] -= 1.0;
|
|---|
| 885 | sRold[d] -= 1.0;
|
|---|
| 886 | sRoldold[d] -= 1.0;
|
|---|
| 887 | }
|
|---|
| 888 | }
|
|---|
| 889 | RMat33Vec3(R, P->Lat.RealBasis, sR);
|
|---|
| 890 | RMat33Vec3(R_old, P->Lat.RealBasis, sRold);
|
|---|
| 891 | RMat33Vec3(R_old_old, P->Lat.RealBasis, sRoldold);
|
|---|
| [6590194] | 892 | //}
|
|---|
| [a0bcf1] | 893 | }
|
|---|
| 894 | }
|
|---|
| 895 | }
|
|---|
| 896 |
|
|---|
| [13d7f6] | 897 | /** Resets the forces array.
|
|---|
| 898 | * \param *P Problem at hand
|
|---|
| 899 | */
|
|---|
| 900 | void ResetForces(struct Problem *P)
|
|---|
| 901 | {
|
|---|
| 902 | struct Ions *I = &P->Ion;
|
|---|
| 903 | double *FIon, *FIonL, *FIonNL, *FEwald;
|
|---|
| 904 | int is, ia, d;
|
|---|
| 905 | for (is=0; is < I->Max_Types; is++) { // go through all types ...
|
|---|
| 906 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) { // .. and each ion of the type
|
|---|
| 907 | FIon = &I->I[is].FIon[NDIM*ia];
|
|---|
| 908 | FIonL = &I->I[is].FIonL[NDIM*ia];
|
|---|
| 909 | FIonNL = &I->I[is].FIonNL[NDIM*ia];
|
|---|
| 910 | FEwald = &I->I[is].FEwald[NDIM*ia];
|
|---|
| 911 | for (d=0; d<NDIM; d++) {
|
|---|
| 912 | FIon[d] = 0.; // zero all as a sign that's been moved
|
|---|
| 913 | FIonL[d] = 0.;
|
|---|
| 914 | FIonNL[d] = 0.;
|
|---|
| 915 | FEwald[d] = 0.;
|
|---|
| 916 | }
|
|---|
| 917 | }
|
|---|
| 918 | }
|
|---|
| 919 | };
|
|---|
| 920 |
|
|---|
| [a0bcf1] | 921 | /** Changes ion's velocity according to acting force.
|
|---|
| 922 | * IonType::U is changed by the weighted force of actual step and last one
|
|---|
| 923 | * according to Newton.
|
|---|
| 924 | * \param *P Problem at hand
|
|---|
| 925 | * \sa UpdateIonsR
|
|---|
| 926 | * \sa CalculateIonForce
|
|---|
| 927 | * \warning R is not changed here.
|
|---|
| 928 | */
|
|---|
| 929 | void UpdateIonsU(struct Problem *P)
|
|---|
| 930 | {
|
|---|
| 931 | struct Ions *I = &P->Ion;
|
|---|
| 932 | int is, ia, d;
|
|---|
| 933 | double *FIon, *FIon_old, *U;
|
|---|
| 934 | double IonMass, a;
|
|---|
| 935 | const int delta_t = P->R.delta_t;
|
|---|
| 936 | for (is=0; is < I->Max_Types; is++) {
|
|---|
| 937 | IonMass = I->I[is].IonMass;
|
|---|
| 938 | if (IonMass < MYEPSILON) fprintf(stderr,"UpdateIonsU: IonMass = %lg",IonMass);
|
|---|
| 939 | a = delta_t*0.5/IonMass; // (F+F_old)/2m = a and thus: v = (F+F_old)/2m * t = (F + F_old) * a
|
|---|
| 940 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
|---|
| 941 | U = &I->I[is].U[NDIM*ia];
|
|---|
| [05fec0] | 942 | if (I->I[is].Constraints[ia] != NULL) { // override current resulting R and U if constrained is given
|
|---|
| 943 | fprintf(stderr, "(%i) Setting constrained motion U at step %d on ion (%i,%i)\n", P->Par.me, I->I[is].Constraints[ia]->step, is, ia);
|
|---|
| 944 | if ((I->I[is].Constraints[ia]->step != P->R.OuterStep) && (P->Par.me == 0))
|
|---|
| 945 | fprintf(stderr, "(%i) WARNING: constraint step %d does not match Outerstep %d for ion (%d,%d)\n", P->Par.me, I->I[is].Constraints[ia]->step, P->R.OuterStep, is, ia);
|
|---|
| 946 | for (d=0; d<NDIM; d++)
|
|---|
| 947 | U[d] = I->I[is].Constraints[ia]->U[d];
|
|---|
| 948 | RemoveConstraintItem(&I->I[is],ia); // and remove this first item from the list
|
|---|
| 949 | } else {
|
|---|
| 950 | FIon = &I->I[is].FIon[NDIM*ia];
|
|---|
| 951 | FIon_old = &I->I[is].FIon_old[NDIM*ia];
|
|---|
| 952 | //if (I->I[is].IMT[ia] == MoveIon)
|
|---|
| 953 | for (d=0; d<NDIM; d++) {
|
|---|
| 954 | U[d] += a * (FIon[d]+FIon_old[d]);
|
|---|
| 955 | }
|
|---|
| 956 | }
|
|---|
| [a0bcf1] | 957 | }
|
|---|
| 958 | }
|
|---|
| 959 | }
|
|---|
| 960 |
|
|---|
| 961 | /** CG process to optimise structure.
|
|---|
| 962 | * \param *P Problem at hand
|
|---|
| 963 | */
|
|---|
| 964 | void UpdateIons(struct Problem *P)
|
|---|
| 965 | {
|
|---|
| 966 | struct Ions *I = &P->Ion;
|
|---|
| 967 | int is, ia, d;
|
|---|
| 968 | double *R, *R_old, *R_old_old, *FIon, *S, *GammaA;
|
|---|
| 969 | double IonFac, GammaB; //, GammaT;
|
|---|
| 970 | double FNorm, StepNorm;
|
|---|
| 971 | for (is=0; is < I->Max_Types; is++) {
|
|---|
| 972 | IonFac = I->I[is].IonFac;
|
|---|
| 973 | GammaA = I->I[is].GammaA;
|
|---|
| 974 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
|---|
| 975 | FIon = &I->I[is].FIon[NDIM*ia];
|
|---|
| 976 | S = &I->I[is].SearchDir[NDIM*ia];
|
|---|
| 977 | GammaB = GammaA[ia];
|
|---|
| 978 | GammaA[ia] = RSP3(FIon,S);
|
|---|
| 979 | FNorm = sqrt(RSP3(FIon,FIon));
|
|---|
| 980 | StepNorm = log(1.+IonFac*FNorm); // Fix Hack
|
|---|
| 981 | if (FNorm != 0)
|
|---|
| 982 | StepNorm /= sqrt(RSP3(FIon,FIon));
|
|---|
| 983 | else
|
|---|
| 984 | StepNorm = 0;
|
|---|
| 985 | //if (GammaB != 0.0) {
|
|---|
| 986 | // GammaT = GammaA[ia]/GammaB;
|
|---|
| 987 | // for (d=0; d>NDIM; d++)
|
|---|
| 988 | // S[d] = FIon[d]+GammaT*S[d];
|
|---|
| 989 | // } else {
|
|---|
| 990 | for (d=0; d<NDIM; d++)
|
|---|
| 991 | S[d] = StepNorm*FIon[d];
|
|---|
| 992 | // }
|
|---|
| 993 | R = &I->I[is].R[NDIM*ia];
|
|---|
| 994 | R_old = &I->I[is].R_old[NDIM*ia];
|
|---|
| 995 | R_old_old = &I->I[is].R_old_old[NDIM*ia];
|
|---|
| 996 | if (I->I[is].IMT[ia] == MoveIon)
|
|---|
| 997 | for (d=0; d<NDIM; d++) {
|
|---|
| 998 | R_old_old[d] = R_old[d];
|
|---|
| 999 | R_old[d] = R[d];
|
|---|
| 1000 | R[d] += S[d]; // FixHack*IonFac;
|
|---|
| 1001 | }
|
|---|
| 1002 | }
|
|---|
| 1003 | }
|
|---|
| 1004 | }
|
|---|
| 1005 |
|
|---|
| 1006 | /** Print coordinates of all ions to stdout.
|
|---|
| 1007 | * \param *P Problem at hand
|
|---|
| [774ae8] | 1008 | * \param actual (1 - don't at current RunStruct#OuterStep, 0 - do)
|
|---|
| [a0bcf1] | 1009 | */
|
|---|
| [774ae8] | 1010 | void OutputIonCoordinates(struct Problem *P, int actual)
|
|---|
| [a0bcf1] | 1011 | {
|
|---|
| [774ae8] | 1012 | //struct RunStruct *R = &P->R;
|
|---|
| [a0bcf1] | 1013 | struct Ions *I = &P->Ion;
|
|---|
| [c510a7] | 1014 | char filename[MAXSTRINGSIZE];
|
|---|
| [774ae8] | 1015 | FILE *output;
|
|---|
| 1016 | int is, ia, nr = 0;
|
|---|
| 1017 | double Bohr = (I->IsAngstroem) ? 1./ANGSTROEMINBOHRRADIUS : 1.;
|
|---|
| 1018 | /* if (P->Par.me == 0 && P->Call.out[ReadOut]) {
|
|---|
| 1019 | fprintf(stderr, "(%i) ======= Differential Ion Coordinates =========\n",P->Par.me);
|
|---|
| [a0bcf1] | 1020 | for (is=0; is < I->Max_Types; is++)
|
|---|
| 1021 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++)
|
|---|
| 1022 | //fprintf(stderr, "(%i) R[%i/%i][%i/%i] = (%e,%e,%e)\n", P->Par.me, is, I->Max_Types, ia, I->I[is].Max_IonsOfType, I->I[is].R[NDIM*ia+0],I->I[is].R[NDIM*ia+1],I->I[is].R[NDIM*ia+2]);
|
|---|
| [774ae8] | 1023 | fprintf(stderr, "Ion_Type%i_%i\t%.6f\t%.6f\t%.6f\t0\t# Atom from StructOpt\n", is+1, ia+1, I->I[is].R[NDIM*ia+0]-I->I[is].R_old[NDIM*ia+0],I->I[is].R[NDIM*ia+1]-I->I[is].R_old[NDIM*ia+1],I->I[is].R[NDIM*ia+2]-I->I[is].R_old[NDIM*ia+2]);
|
|---|
| [a0bcf1] | 1024 | fprintf(stderr, "(%i) =========================================\n",P->Par.me);
|
|---|
| [774ae8] | 1025 | }*/
|
|---|
| 1026 | if (actual)
|
|---|
| 1027 | sprintf(filename, "%s.MD", P->Call.MainParameterFile);
|
|---|
| 1028 | else
|
|---|
| 1029 | sprintf(filename, "%s.opt", P->Call.MainParameterFile);
|
|---|
| 1030 | if (((actual) && (P->R.OuterStep <= 0)) || ((!actual) && (P->R.StructOptStep == 1))) { // on first step write complete config file
|
|---|
| 1031 | WriteParameters(P, filename);
|
|---|
| 1032 | } else { // afterwards simply add the current ion coordinates as constrained steps
|
|---|
| 1033 | if ((P->Par.me == 0) && (output = fopen(filename,"a"))) {
|
|---|
| 1034 | fprintf(output,"# ====== MD step %d ========= \n", (actual) ? P->R.OuterStep : P->R.StructOptStep);
|
|---|
| 1035 | for (is=0; is < I->Max_Types; is++)
|
|---|
| 1036 | for (ia=0;ia<I->I[is].Max_IonsOfType;ia++) {
|
|---|
| 1037 | fprintf(output,"Ion_Type%i_%i\t%2.9f\t%2.9f\t%2.9f\t%i\t%2.9f\t%2.9f\t%2.9f\t# Number in molecule %i\n", is+1, ia+1, I->I[is].R[0+NDIM*ia]*Bohr, I->I[is].R[1+NDIM*ia]*Bohr, I->I[is].R[2+NDIM*ia]*Bohr, I->I[is].IMT[ia], I->I[is].U[0+NDIM*ia]*Bohr, I->I[is].U[1+NDIM*ia]*Bohr, I->I[is].U[2+NDIM*ia]*Bohr, ++nr);
|
|---|
| 1038 | }
|
|---|
| 1039 | fflush(output);
|
|---|
| 1040 | }
|
|---|
| [a0bcf1] | 1041 | }
|
|---|
| 1042 | }
|
|---|
| 1043 |
|
|---|
| 1044 | /** Calculates kinetic energy (Ions::EKin) of movable Ions.
|
|---|
| 1045 | * Does 0.5 / IonType::IonMass * IonTye::U^2 for each Ion,
|
|---|
| 1046 | * also retrieves actual temperatur by 2/3 from just
|
|---|
| 1047 | * calculated Ions::EKin.
|
|---|
| 1048 | * \param *P Problem at hand
|
|---|
| 1049 | */
|
|---|
| 1050 | void CalculateEnergyIonsU(struct Problem *P)
|
|---|
| 1051 | {
|
|---|
| 1052 | struct Ions *I = &P->Ion;
|
|---|
| 1053 | int is, ia, d;
|
|---|
| 1054 | double *U;
|
|---|
| 1055 | double IonMass, a, ekin = 0;
|
|---|
| 1056 | for (is=0; is < I->Max_Types; is++) {
|
|---|
| 1057 | IonMass = I->I[is].IonMass;
|
|---|
| 1058 | a = 0.5*IonMass;
|
|---|
| 1059 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
|---|
| 1060 | U = &I->I[is].U[NDIM*ia];
|
|---|
| [3716b2] | 1061 | //if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
|
|---|
| [a0bcf1] | 1062 | for (d=0; d<NDIM; d++) {
|
|---|
| 1063 | ekin += a * U[d]*U[d];
|
|---|
| 1064 | }
|
|---|
| 1065 | }
|
|---|
| 1066 | }
|
|---|
| 1067 | I->EKin = ekin;
|
|---|
| 1068 | I->ActualTemp = (2./(3.*I->Max_TotalIons)*I->EKin);
|
|---|
| 1069 | }
|
|---|
| 1070 |
|
|---|
| 1071 | /** Scales ion velocities to match temperature.
|
|---|
| 1072 | * In order to match Ions::ActualTemp with desired RunStruct::TargetTemp
|
|---|
| 1073 | * each velocity in each dimension (for all types, all ions) is scaled
|
|---|
| 1074 | * with the ratio of these two. Ions::EKin is recalculated.
|
|---|
| 1075 | * \param *P Problem at hand
|
|---|
| 1076 | */
|
|---|
| 1077 | void ScaleTemp(struct Problem *P)
|
|---|
| 1078 | {
|
|---|
| 1079 | struct Ions *I = &P->Ion;
|
|---|
| 1080 | int is, ia, d;
|
|---|
| 1081 | double *U;
|
|---|
| 1082 | double IonMass, a, ekin = 0;
|
|---|
| 1083 | if (I->ActualTemp < MYEPSILON) fprintf(stderr,"ScaleTemp: I->ActualTemp = %lg",I->ActualTemp);
|
|---|
| 1084 | double ScaleTempFactor = sqrt(P->R.TargetTemp/I->ActualTemp);
|
|---|
| 1085 | for (is=0; is < I->Max_Types; is++) {
|
|---|
| 1086 | IonMass = I->I[is].IonMass;
|
|---|
| 1087 | a = 0.5*IonMass;
|
|---|
| 1088 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
|---|
| 1089 | U = &I->I[is].U[NDIM*ia];
|
|---|
| [3716b2] | 1090 | //if (I->I[is].IMT[ia] == MoveIon) // even FixedIon moves, only not by other's forces
|
|---|
| [a0bcf1] | 1091 | for (d=0; d<NDIM; d++) {
|
|---|
| 1092 | U[d] *= ScaleTempFactor;
|
|---|
| 1093 | ekin += a * U[d]*U[d];
|
|---|
| 1094 | }
|
|---|
| 1095 | }
|
|---|
| 1096 | }
|
|---|
| 1097 | I->EKin = ekin;
|
|---|
| 1098 | I->ActualTemp = (2./(3.*I->Max_TotalIons)*I->EKin);
|
|---|
| 1099 | }
|
|---|
| 1100 |
|
|---|
| 1101 | /** Calculates mean force vector as stop criteria in structure optimization.
|
|---|
| 1102 | * Calculates a mean force vector per ion as the usual euclidian norm over total
|
|---|
| 1103 | * number of ions and dimensions, being the sum of each ion force (for all type,
|
|---|
| 1104 | * all ions, all dims) squared.
|
|---|
| 1105 | * The mean force is archived in RunStruct::MeanForce and printed to screen.
|
|---|
| 1106 | * \param *P Problem at hand
|
|---|
| 1107 | */
|
|---|
| 1108 | void GetOuterStop(struct Problem *P)
|
|---|
| 1109 | {
|
|---|
| 1110 | struct RunStruct *R = &P->R;
|
|---|
| 1111 | struct Ions *I = &P->Ion;
|
|---|
| 1112 | int is, ia, IonNo=0, i, d;
|
|---|
| 1113 | double MeanForce = 0.0;
|
|---|
| 1114 | for (is=0; is < I->Max_Types; is++)
|
|---|
| 1115 | for (ia=0; ia < I->I[is].Max_IonsOfType; ia++) {
|
|---|
| 1116 | IonNo++;
|
|---|
| 1117 | for (d=0; d <NDIM; d++)
|
|---|
| 1118 | MeanForce += I->I[is].FIon[d+NDIM*ia]*I->I[is].FIon[d+NDIM*ia];
|
|---|
| 1119 | }
|
|---|
| 1120 | for (i=MAXOLD-1; i > 0; i--)
|
|---|
| 1121 | R->MeanForce[i] = R->MeanForce[i-1];
|
|---|
| 1122 | MeanForce = sqrt(MeanForce/(IonNo*NDIM));
|
|---|
| 1123 | R->MeanForce[0] = MeanForce;
|
|---|
| 1124 | if (P->Call.out[LeaderOut] && (P->Par.me == 0))
|
|---|
| 1125 | fprintf(stderr,"MeanForce = %e\n", R->MeanForce[0]);
|
|---|
| 1126 | }
|
|---|
| 1127 |
|
|---|
| 1128 |
|
|---|