Octopus
io_function.F90
Go to the documentation of this file.
1!! Copyright (C) 2002-2006 M. Marques, A. Castro, A. Rubio, G. Bertsch
2!!
3!! This program is free software; you can redistribute it and/or modify
4!! it under the terms of the GNU General Public License as published by
5!! the Free Software Foundation; either version 2, or (at your option)
6!! any later version.
7!!
8!! This program is distributed in the hope that it will be useful,
9!! but WITHOUT ANY WARRANTY; without even the implied warranty of
10!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
11!! GNU General Public License for more details.
12!!
13!! You should have received a copy of the GNU General Public License
14!! along with this program; if not, write to the Free Software
15!! Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
16!! 02110-1301, USA.
17!!
18
19#include "global.h"
20
23 use atom_oct_m
24 use box_oct_m
26 use comm_oct_m
28 use cube_oct_m
29 use debug_oct_m
31 use fft_oct_m
32 use global_oct_m
33 use index_oct_m
34 use io_oct_m
36 use io_csv_oct_m
37 use iso_c_binding
40 use mesh_oct_m
43 use mpi_oct_m
46#if defined(HAVE_NETCDF)
47 use netcdf
48#endif
50 use parser_oct_m
53 use space_oct_m
55 use unit_oct_m
57 use utils_oct_m
59 use vtk_oct_m
60
61 implicit none
62
63 private
64 public :: &
80
81#if defined(HAVE_NETCDF)
82 public :: &
85#endif
86
88 integer, parameter, private :: &
89 DOUTPUT_KIND = 1, &
90 zoutput_kind = -1
91
93 character(len=3), parameter :: &
94 index2label(3) = (/ 're ', 'im ', 'abs' /)
95
96
99 end interface io_function_output_vector
100
103 end interface io_function_output_vector_bz
104
107 end interface io_function_output_global_bz
108
111 end interface io_function_output_supercell
112
113
115contains
116
117 ! ---------------------------------------------------------
118 subroutine io_function_read_what_how_when(namespace, space, what, how, output_interval, &
119 what_tag_in, how_tag_in, output_interval_tag_in, ignore_error)
120 type(namespace_t), intent(in) :: namespace
121 class(space_t), intent(in) :: space
122 logical, intent(inout) :: what(MAX_OUTPUT_TYPES)
123 integer(int64), intent(out) :: how(0:MAX_OUTPUT_TYPES)
124 integer, intent(out) :: output_interval(0:MAX_OUTPUT_TYPES)
125 character(len=*), optional, intent(in) :: what_tag_in
126 character(len=*), optional, intent(in) :: how_tag_in
127 character(len=*), optional, intent(in) :: output_interval_tag_in
128 logical, optional, intent(in) :: ignore_error
129
130 type(block_t) :: blk
131 integer :: ncols, nrows, iout, column_index, what_i
132 integer(int64) :: what_no_how(12)
133 character(len=80) :: what_tag
134 character(len=80) :: how_tag
135 character(len=80) :: output_interval_tag
136 character(len=80) :: output_column_marker
138
139 how = 0_8
140 output_interval = 0_8
141 ncols = 0
142 what_tag = optional_default(what_tag_in, 'Output')
143 how_tag = optional_default(how_tag_in, 'OutputFormat')
144 output_interval_tag = optional_default(output_interval_tag_in, 'OutputInterval')
145
146 call messages_obsolete_variable(namespace, 'OutputHow', 'OutputFormat')
147 call messages_obsolete_variable(namespace, 'Output_KPT', 'Output')
148 call messages_obsolete_variable(namespace, 'OutputLDA_U' , 'Output')
149 call messages_obsolete_variable(namespace, 'OutputEvery', 'OutputInterval/RestartWriteInterval')
150
151 !%Variable Output
152 !%Type block
153 !%Default none
154 !%Section Output
155 !%Description
156 !% Specifies what to print.
157 !% Each output must be in a separate row. Optionally individual output formats and output intervals can be defined
158 !% for each row or they can be read separately from <tt>OutputFormat</tt> and <tt>OutputInterval</tt> variables
159 !% in the input file.
160 !% The output files are written at the end of the run into the output directory for the
161 !% relevant kind of run (<i>e.g.</i> <tt>static</tt> for <tt>CalculationMode = gs</tt>).
162 !% Time-dependent simulations print only per iteration, including always the last. The frequency of output per iteration
163 !% (available for <tt>CalculationMode</tt> = <tt>gs</tt>, <tt>unocc</tt>, <tt>td</tt>, and <tt>opt_control</tt>)
164 !% is set by <tt>OutputInterval</tt> and the directory is set by <tt>OutputIterDir</tt>.
165 !% For linear-response run modes, the derivatives of many quantities can be printed, as listed in
166 !% the options below. Indices in the filename are labelled as follows:
167 !% <tt>sp</tt> = spin (or spinor component), <tt>k</tt> = <i>k</i>-point, <tt>st</tt> = state/band.
168 !% There is no tag for directions, given as a letter. The perturbation direction is always
169 !% the last direction for linear-response quantities, and a following +/- indicates the sign of the frequency.
170 !%
171 !% Example (minimal):
172 !% <br><br><tt>%Output
173 !% <br>&nbsp;&nbsp;density
174 !% <br>&nbsp;&nbsp;potential
175 !% <br>%<br></tt>
176 !%
177 !% Example (with OutputFormat):
178 !% <br><br><tt>%Output
179 !% <br>&nbsp;&nbsp;density | cube + axis_z
180 !% <br>&nbsp;&nbsp;potential | cube
181 !% <br>%<br></tt>
182 !%
183 !% Example (with OutputFormat, incomplete):
184 !% <br><br><tt>%Output
185 !% <br>&nbsp;&nbsp;density | cube + axis_z
186 !% <br>&nbsp;&nbsp;potential
187 !% <br>%<br></tt>
188 !%
189 !% Example (tagged):
190 !% <br><br><tt>%Output
191 !% <br>&nbsp;&nbsp;density | "output_format" | cube + axis_z | "output_interval" | 50
192 !% <br>&nbsp;&nbsp;potential | "output_format" | cube | "output_interval" | 20
193 !% <br>%<br></tt>
194 !%
195 !% Example (tagged, incomplete):
196 !% <br><br><tt>%Output
197 !% <br>&nbsp;&nbsp;density | "output_format" | cube + axis_z
198 !% <br>&nbsp;&nbsp;potential | "output_interval" | 20
199 !% <br>%<br></tt>
200 !% Missing information for the incomplete blocks will be parsed form the out-of-block
201 !% definitions. It is also possible to mix the order of columns in the tagged format.
202 !% See <tt>OutputFormat</tt>, and <tt>OutputInterval</tt>.
203 !%Option potential 1
204 !% Outputs Kohn-Sham potential, separated by parts. File names are <tt>v0</tt> for
205 !% the local part of the ionic potential, <tt>vc</tt> for the classical potential (if it exists),
206 !% <tt>vh</tt> for the Hartree potential, <tt>vks</tt> for the local part of the Kohn-Sham potential, and
207 !% <tt>vxc-</tt> for the exchange-correlation potentials. For <tt>vks</tt> and <tt>vxc</tt>,
208 !% a suffix for spin is added in the spin-polarized case.
209 !%Option density 2
210 !% Outputs density. The output file is called <tt>density-</tt>, or <tt>lr_density-</tt> in linear response.
211 !%Option wfs 3
212 !% Outputs wavefunctions. Which wavefunctions are to be printed is specified
213 !% by the variable <tt>OutputWfsNumber</tt> -- see below. The output file is called
214 !% <tt>wf-</tt>, or <tt>lr_wf-</tt> in linear response.
215 !%Option wfs_sqmod 4
216 !% Outputs modulus squared of the wavefunctions.
217 !% The output file is called <tt>sqm-wf-</tt>. For linear response, the filename is <tt>sqm_lr_wf-</tt>.
218 !%Option geometry 5
219 !% Outputs file containing the coordinates of the atoms treated within quantum mechanics.
220 !% If <tt>OutputFormat = xyz</tt>, the file is called <tt>geometry.xyz</tt>; a
221 !% file <tt>crystal.xyz</tt> is written with a supercell geometry if the system is periodic;
222 !% if point charges were defined in the PDB file (see <tt>PDBCoordinates</tt>), they will be output
223 !% in the file <tt>geometry_classical.xyz</tt>.
224 !% If <tt>OutputFormat = xcrysden</tt>, a file called <tt>geometry.xsf</tt> is written.
225 !%Option current 6
226 !% Outputs the total current density. The output file is called <tt>current-</tt>.
227 !% For linear response, the filename is <tt>lr_current-</tt>.
228 !%Option ELF 7
229 !% Outputs electron localization function (ELF). The output file is called <tt>elf-</tt>,
230 !% or <tt>lr_elf-</tt> in linear response, in which case the associated function D is also written,
231 !% as <tt>lr_elf_D-</tt>. Only in 2D and 3D.
232 !%Option ELF_basins 8
233 !% Outputs basins of attraction of the ELF. The output file is called
234 !% <tt>elf_rs_basins.info</tt>. Only in 2D and 3D.
235 !%Option Bader 9
236 !% Outputs Laplacian of the density which shows lone pairs, bonded charge concentrations
237 !% and regions subject to electrophilic or nucleophilic attack.
238 !% See RF Bader, <i>Atoms in Molecules: A Quantum Theory</i> (Oxford Univ. Press, Oxford, 1990).
239 !%Option el_pressure 10
240 !% Outputs electronic pressure. See Tao, Vignale, and Tokatly, <i>Phys Rev Lett</i> <b>100</b>, 206405 (2008).
241 !%Option matrix_elements 11
242 !% Outputs a series of matrix elements of the Kohn-Sham states. What is output can
243 !% be controlled by the <tt>OutputMatrixElements</tt> variable.
244 !%Option pol_density 12
245 !% Outputs dipole-moment density <tt>dipole_density-</tt>, or polarizability density <tt>alpha_density-</tt>
246 !% in linear response. If <tt>ResponseMethod = finite_differences</tt>, the hyperpolarizability density
247 !% <tt>beta_density-</tt> is also printed.
248 !%Option mesh_r 13
249 !% Outputs values of the coordinates over the grid. Files
250 !% will be called <tt>mesh_r-</tt> followed by the direction.
251 !%Option kinetic_energy_density 14
252 !% Outputs kinetic-energy density, defined as:
253 !%
254 !% <math>\tau_\sigma(\vec{r}) = \sum_{i=1}^{N_\sigma}
255 !% \left| \vec{\nabla} \phi_{i\sigma}(\vec{r}) \right|^2\,. </math>
256 !%
257 !% The index <math>\sigma</math> is the spin index for the spin-polarized case,
258 !% or if you are using spinors. For spin-unpolarized calculations, you
259 !% get the total kinetic-energy density. The previous expression assumes full
260 !% or null occupations. If fractional occupation numbers, each term in the sum
261 !% is weighted by the occupation. Also, if we are working with an infinite
262 !% system, all <i>k</i>-points are summed up, with their corresponding weights. The
263 !% files will be called <tt>tau-sp1</tt> and <tt>tau-sp2</tt>, if the spin-resolved kinetic
264 !% energy density is produced (runs in spin-polarized and spinors mode), or
265 !% only <tt>tau</tt> if the run is in spin-unpolarized mode.
266 !%Option dos 15
267 !% Outputs density of states. See <tt>DOSEnergyMax</tt>, <tt>DOSEnergyMin</tt>, <tt>DOSEnergyPoints</tt>,
268 !% and <tt>DOSGamma</tt>.
269 !%Option tpa 16
270 !% Outputs transition-potential approximation (TPA) matrix elements, using <math>\vec{q}</math>-vector specified
271 !% by <tt>MomentumTransfer</tt>.
272 !%Option forces 17
273 !% Outputs file <tt>forces.xsf</tt> containing structure and forces on the atoms as
274 !% a vector associated with each atom, which can be visualized with XCrySDen.
275 !%Option wfs_fourier 18
276 !% (Experimental) Outputs wavefunctions in Fourier space. This is
277 !% only implemented for the ETSF file format output. The file will
278 !% be called <tt>wfs-pw-etsf.nc</tt>.
279 !%Option xc_density 19
280 !% Outputs the XC density, which is the charge density that
281 !% generates the XC potential. (This is <math>-1/4\pi</math> times
282 !% the Laplacian of the XC potential). The files are called <tt>nxc</tt>.
283 !%Option PES_wfs 20
284 !% Outputs the photoelectron wavefunctions. The file name is <tt>pes_wfs-</tt>
285 !% plus the orbital number.
286 !%Option PES_density 21
287 !% Outputs the photolectron density. Output file is <tt>pes_dens-</tt> plus spin species if
288 !% spin-polarized calculation is performed.
289 !%Option PES 22
290 !% Outputs the time-dependent photoelectron spectrum.
291 !%Option BerkeleyGW 23
292 !% Output for a run with <a href=http://www.berkeleygw.org>BerkeleyGW</a>.
293 !% See <tt>Output::BerkeleyGW</tt> for further specification.
294 !%Option delta_perturbation 24
295 !% Outputs the "kick", or time-delta perturbation applied to compute optical response in real time.
296 !%Option external_td_potential 25
297 !% Outputs the (scalar) time-dependent potential.
298 !%Option mmb_wfs 26
299 !% Triggers the ModelMB wavefunctions to be output for each state.
300 !%Option mmb_den 27
301 !% Triggers the ModelMB density matrix to be output for each state, and the particles
302 !% specified by the <tt>DensitytoCalc</tt> block. Calculates, and outputs, the reduced density
303 !% matrix. For the moment the trace is made over the second dimension, and
304 !% the code is limited to 2D. The idea is to model <i>N</i> particles in 1D as an
305 !% <i>N</i>-dimensional non-interacting problem, then to trace out <i>N</i>-1 coordinates.
306 !%Option potential_gradient 28
307 !% Prints the gradient of the potential.
308 !%Option energy_density 29
309 !% Outputs the total energy density to a file called
310 !% <tt>energy_density</tt>.
311 !%Option heat_current 30
312 !% Outputs the total heat current density. The output file is
313 !% called <tt>heat_current-</tt>.
314 !%Option photon_correlator 31
315 !% Outputs the electron-photon correlation function. The output file is
316 !% called <tt>photon_correlator</tt>.
317 !%Option J_flow 32
318 !% todo: document J_flow option!
319 !%Option current_kpt 33
320 !% Outputs the current density resolved in momentum space. The output file is called <tt>current_kpt-</tt>.
321 !%Option density_kpt 34
322 !% Outputs the electronic density resolved in momentum space.
323 !%Option occ_matrices 35
324 !% Only for DFT+U calculations.
325 !% Outputs the occupation matrices of DFT+U
326 !%Option effectiveU 36
327 !% Only for DFT+U calculations.
328 !% Outputs the value of the effectiveU for each atoms
329 !%Option magnetization 37
330 !% Only for DFT+U calculations.
331 !% Outputs file containing structure and magnetization of the localized subspace
332 !% on the atoms as a vector associated with each atom, which can be visualized.
333 !% For the moment, it only works if a +U is added on one type of orbital per atom.
334 !%Option local_orbitals 38
335 !% Only for DFT+U calculations.
336 !% Outputs the localized orbitals that form the correlated subspace
337 !%Option kanamoriU 39
338 !% Only for DFT+U calculations.
339 !% Outputs the Kanamori interaction parameters U, U`, and J.
340 !% These parameters are not determined self-consistently, but are taken from the
341 !% occupation matrices and Coulomb integrals comming from a standard +U calculation.
342 !%Option xc_torque 40
343 !% Outputs the exchange-correlation torque. Only for the spinor case and in the 3D case.
344 !%Option eigenval_kpt 41
345 !% Outputs the eigenvalues resolved in momentum space, with one file for each band.
346 !%Option stress 42
347 !% Outputs the stress tensor and each of its contributing terms
348 !%Option current_dia 43
349 !% Outputs the diamagnetic current density from a non-uniform vector potential.
350 !% The output file is called <tt>current_dia-</tt>.
351 !%Option jdos 44
352 !% Outputs the joint density of states.
353 !% The same variables as for the regular DOS are used to control the energies and broadening.
354 !%Option ldos 45
355 !% Outputs the local density of states
356 !% The broadening uses the same value as the regular DOS. The energies at which the LDOS is computed are
357 !% specified by the <tt>LDOSEnergies</tt> block.
358 !%End
359
360 !%Variable OutputFormat
361 !%Type flag
362 !%Default 0
363 !%Section Output
364 !%Description
365 !% Describes the format of the output files.
366 !% This variable can also be defined inside the <tt>Output</tt> block.
367 !% See <tt>Output</tt>.
368 !% Example: <tt>axis_x + plane_x + dx</tt>
369 !%Option axis_x bit(0)
370 !% The values of the function on the <i>x</i> axis are printed. The string <tt>.y=0,z=0</tt> is appended
371 !% to previous file names.
372 !%Option axis_y bit(1)
373 !% The values of the function on the <i>y</i> axis are printed. The string <tt>.x=0,z=0</tt> is appended
374 !% to previous file names.
375 !%Option axis_z bit(2)
376 !% The values of the function on the <i>z</i> axis are printed. The string <tt>.x=0,y=0</tt> is appended
377 !% to previous file names.
378 !%Option plane_x bit(3)
379 !% A plane slice at <i>x</i> = 0 is printed. The string <tt>.x=0</tt> is appended
380 !% to previous file names.
381 !%Option plane_y bit(4)
382 !% A plane slice at <i>y</i> = 0 is printed. The string <tt>.y=0</tt> is appended
383 !% to previous file names.
384 !%Option plane_z bit(5)
385 !% A plane slice at <i>z</i> = 0 is printed. The string <tt>.z=0</tt> is appended to
386 !% previous file names.
387 !%Option dx bit(6)
388 !% For printing three-dimensional information, the open-source program
389 !% visualization tool <a href=http://www.opendx.org>OpenDX</a> can be used. The string
390 !% <tt>.dx</tt> is appended to previous file names. Available only in 3D.
391 !%Option netcdf bit(7)
392 !% Outputs in <a href=http://www.unidata.ucar.edu/packages/netcdf>NetCDF</a> format. This file
393 !% can then be read, for example, by OpenDX. The string <tt>.ncdf</tt> is appended to previous file names.
394 !% Requires the NetCDF library. Only writes the real part of complex functions.
395 !%Option mesh_index bit(8)
396 !% Generates output files of a given quantity (density, wavefunctions, ...) which include
397 !% the internal numbering of mesh points. Since this mode produces large datafiles this is only
398 !% useful for small meshes and debugging purposes.
399 !% The output can also be used to display the mesh directly. A Gnuplot script for mesh visualization
400 !% can be found under <tt>PREFIX/share/octopus/util/display_mesh_index.gp</tt>.
401 !%Option xcrysden bit(9)
402 !% A format for printing structures and three-dimensional information, which can be visualized by
403 !% the free open-source program <a href=http://www.xcrysden.org>XCrySDen</a> and others. The string
404 !% <tt>.xsf</tt> is appended to previous file names. Note that lattice vectors and coordinates are as
405 !% specified by <tt>UnitsOutput</tt>. Available in 2D and 3D.
406 !%Option matlab bit(10)
407 !% In combination with <tt>plane_x</tt>, <tt>plane_y</tt> and
408 !% <tt>plane_z</tt>, this option produces output files which are
409 !% suitable for 2D Matlab functions like <tt>mesh()</tt>,
410 !% <tt>surf()</tt>, or <tt>waterfall()</tt>. To load these files
411 !% into Matlab you can use, <i>e.g.</i>
412 !%<tt>
413 !% >> density = load('static/density-1.x=0.matlab.abs');
414 !% >> mesh(density);
415 !%</tt>
416 !%Option meshgrid bit(11)
417 !% Outputs in Matlab mode the internal mesh in a format similar to
418 !%<tt>
419 !% >> [x,y] = meshgrid(-2:.2:2,-1:.15:1)
420 !%</tt>
421 !% The <i>x</i> meshgrid is contained in a file <tt>*.meshgrid.x</tt> and the <i>y</i>-grid can be found in
422 !% <tt>*.meshgrid.y</tt>.
423 !%Option boundary_points bit(12)
424 !% This option includes the output of the mesh enlargement. Default is without.
425 !% Supported only by <tt>binary</tt>, <tt>axis</tt>, <tt>plane</tt>, <tt>mesh_index</tt>,
426 !% and <tt>matlab</tt> formats.
427 !% Not all types of <tt>Output</tt> will have this information available. Not supported when parallel in domains.
428 !%Option binary bit(13)
429 !% Plain binary, new format.
430 !%Option etsf bit(14)
431 !% <a href=http://www.etsf.eu/resources/software/standardization_project>ETSF file format</a>.
432 !% Requires the ETSF_IO library. Applies only to <tt>Output = density</tt>, <tt>geometry</tt>,
433 !% <tt>wfs</tt>, and/or <tt>wfs_fourier</tt>.
434 !%Option xyz bit(15)
435 !% Geometry will be output in XYZ format. Does not affect other outputs.
436 !%Option cube bit(16)
437 !% Generates output in the <a href=http://paulbourke.net/dataformats/cube>cube file format</a>.
438 !% Available only in 3D. Only writes the real part of complex functions.
439 !% This output format always uses atomic units.
440 !%Option bild bit(19)
441 !% Generates output in <a href=http://plato.cgl.ucsf.edu/chimera/docs/UsersGuide/bild.html>BILD format</a>.
442 !%Option vtk bit(20)
443 !% Generates output in <a href=http://www.vtk.org/VTK/img/file-formats.pdf>VTK legacy format</a>.
444 !%Option integrate_xy bit(21)
445 !% Integrates the function in the x-y plane and the result on the <i>z</i> axis is printed.
446 !%Option integrate_xz bit(22)
447 !% Integrates the function in the x-z plane and the result on the <i>y</i> axis is printed
448 !%Option integrate_yz bit(23)
449 !% Integrates the function in the y-z plane and the result on the <i>x</i> axis is printed
450 !%Option ascii bit(24)
451 !% Plain text format regardless of dimensionality. For the moment only employed by the oct-phototoelectron_spectrum
452 !% post-processing utility.
453 !%End
454
455 !%Variable OutputInterval
456 !%Type integer
457 !%Default 50
458 !%Section Output
459 !%Description
460 !% The output requested by variable <tt>Output</tt> is written
461 !% to the directory <tt>OutputIterDir</tt>
462 !% when the iteration number is a multiple of the <tt>OutputInterval</tt> variable.
463 !% Subdirectories are named Y.X, where Y is <tt>td</tt>, <tt>scf</tt>, or <tt>unocc</tt>, and
464 !% X is the iteration number. To use the working directory, specify <tt>"."</tt>
465 !% (Output of restart files is instead controlled by <tt>RestartWriteInterval</tt>.)
466 !% Must be >= 0. If it is 0, then no output is written. For <tt>gs</tt> and <tt>unocc</tt>
467 !% calculations, <tt>OutputDuringSCF</tt> must be set too for this output to be produced.
468 !% This variable can also be defined inside the <tt>Output</tt> block.
469 !% See <tt>Output</tt>.
470 !%End
471
472
473 what_no_how = (/ option__output__matrix_elements, option__output__berkeleygw, option__output__dos, &
474 option__output__tpa, option__output__mmb_den, option__output__j_flow, option__output__occ_matrices, &
475 option__output__effectiveu, option__output__magnetization, option__output__kanamoriu, option__output__stress, option__output__jdos /)
476
477 if (parse_block(namespace, what_tag, blk) == 0) then
478 nrows = parse_block_n(blk)
479 do iout = 0, nrows - 1
480 ncols = max(ncols , parse_block_cols(blk, iout))
481 end do
482
483 if (ncols == 1) then
484 !new format, Type 0
485 !%Output
486 ! density
487 ! wfs
488 !%
489 do iout = 1, nrows
490 call parse_block_integer(blk, iout - 1, 0, what_i)
491 if (.not. varinfo_valid_option(what_tag, what_i)) then
492 call messages_input_error(namespace, what_tag)
493 end if
494 if (what_i > 0) then
495 what(what_i) = .true.
496 call parse_variable(namespace, output_interval_tag, 50, output_interval(what_i))
497 if (((what_tag == 'Output') .and. (.not. any(what_no_how == what_i)))&
498 .or. (what_tag /= 'Output')) then
499 call parse_variable(namespace, how_tag, 0, how(what_i))
500 if (.not. varinfo_valid_option(how_tag, how(what_i), is_flag=.true.)) then
501 call messages_input_error(namespace, how_tag)
502 end if
503 end if
504 end if
505 end do
506 else if (ncols == 2) then
507 !new format, Type 1
508 !%Output
509 ! density | cube + axis_z
510 ! wfs | cube
511 !%
512
513 do iout = 1, nrows
514 call parse_block_integer(blk, iout - 1, 0, what_i)
515 if (.not. varinfo_valid_option(what_tag, what_i)) then
516 call messages_input_error(namespace, what_tag)
517 end if
518 if (what_i > 0) then
519 what(what_i) = .true.
520 call parse_variable(namespace, output_interval_tag, 50, output_interval(what_i))
521 if (((what_tag == 'Output') .and. (.not. any(what_no_how == what_i)))&
522 .or. (what_tag /= 'Output')) then
523 call parse_block_integer(blk, iout - 1, 1, how(what_i))
524 if (how(what_i) == 0) call parse_variable(namespace, how_tag, 0, how(what_i))
525 if (.not. varinfo_valid_option(how_tag, how(what_i), is_flag=.true.)) then
526 call messages_input_error(namespace, how_tag)
527 end if
528 end if
529 end if
530 end do
531
532 else
533 !new format, Type 2 (tagged)
534 !%Output
535 ! density | "output_interval" | 10 | "output_format" | cube + axis_z
536 ! wfs | "output_format" | cube | "output_interval" | 50
537 !%
538 !
539 ! OR
540 !
541 !%Output
542 ! density | "output_interval" | 10 | "output_format" | cube + axis_z
543 ! wfs | "output_format" | cube
544 !%
545 do iout = 1, nrows
546 call parse_block_integer(blk, iout - 1, 0, what_i)
547 if (.not. varinfo_valid_option(what_tag, what_i)) then
548 call messages_input_error(namespace, what_tag)
549 end if
550 if (what_i > 0) then
551 what(what_i) = .true.
552 do column_index = 0, 1
553 call parse_block_string(blk, iout - 1, 1 + column_index * 2, output_column_marker)
554 if (output_column_marker == 'output_interval') then
555 call parse_block_integer(blk, iout - 1, 2 + column_index * 2, output_interval(what_i))
556 else if (output_column_marker == 'output_format') then
557 if (((what_tag == 'Output') .and. (.not. any(what_no_how == what_i)))&
558 .or. (what_tag /= 'Output')) then
559 call parse_block_integer(blk, iout - 1, 2 + column_index * 2, how(what_i))
560 if (.not. varinfo_valid_option(how_tag, how(what_i), is_flag=.true.)) then
561 call messages_input_error(namespace, how_tag)
562 end if
563 end if
564 else if (len_trim(output_column_marker) /= 0) then
565 ! Unknown output_column_marker
566 call messages_input_error(namespace, what_tag)
567 else
568 ! no output_column_marker -> full output info is not in this block
569 if (output_interval(what_i) == 0) then
570 call parse_variable(namespace, output_interval_tag, 50, output_interval(what_i))
571 end if
572 if (how(what_i) == 0) then
573 if (((what_tag == 'Output') .and. (.not. any(what_no_how == what_i)))&
574 .or. (what_tag /= 'Output')) then
575 call parse_variable(namespace, how_tag, 0, how(what_i))
576 if (.not. varinfo_valid_option(how_tag, how(what_i), is_flag=.true.)) then
577 call messages_input_error(namespace, how_tag)
578 end if
579 end if
580 end if
581 end if
582 end do
583 end if
584 end do
585 end if
586 call parse_block_end(blk)
587 else
588
589 call messages_variable_is_block(namespace, what_tag)
590
591 ! Output block does not exist but we may have OutputHow/OutputInterval
592 call parse_variable(namespace, how_tag, 0, how(0))
593 call parse_variable(namespace, output_interval_tag, 50, output_interval(0))
594 end if
595
596
597 do what_i = lbound(what, 1), ubound(what, 1)
598 if (what_tag == 'Output') then
599 if (what(what_i) .and. (.not. any(what_no_how == what_i))) then
600 if (.not. varinfo_valid_option(how_tag, how(what_i), is_flag=.true.)) then
601 call messages_input_error(namespace, how_tag)
602 end if
603
604 if (how(what_i) == 0 .and. .not. optional_default(ignore_error, .false.)) then
605 write(message(1), '(a)') 'Must specify output method with variable OutputFormat.'
606 call messages_fatal(1, only_root_writes = .true., namespace=namespace)
607 end if
608
609 ! some modes are not available in some circumstances
610 if (space%dim == 1) then
611 if (bitand(how(what_i), option__outputformat__axis_y) /= 0) then
612 message(1) = "OutputFormat = axis_y not available with Dimensions = 1."
613 call messages_fatal(1, namespace=namespace)
614 end if
615 if (bitand(how(what_i), option__outputformat__plane_z) /= 0) then
616 message(1) = "OutputFormat = plane_z not available with Dimensions = 1."
617 call messages_fatal(1, namespace=namespace)
618 end if
619 if (bitand(how(what_i), option__outputformat__xcrysden) /= 0) then
620 message(1) = "OutputFormat = xcrysden not available with Dimensions = 1."
621 call messages_fatal(1, namespace=namespace)
622 end if
623 end if
624
625 if (space%dim <= 2) then
626 if (bitand(how(what_i), option__outputformat__axis_z) /= 0) then
627 message(1) = "OutputFormat = axis_z not available with Dimensions <= 2."
628 call messages_fatal(1, namespace=namespace)
629 end if
630 if (bitand(how(what_i), option__outputformat__plane_x) /= 0) then
631 message(1) = "OutputFormat = plane_x not available with Dimensions <= 2."
632 call messages_fatal(1, namespace=namespace)
633 end if
634 if (bitand(how(what_i), option__outputformat__plane_y) /= 0) then
635 message(1) = "OutputFormat = plane_y not available with Dimensions <= 2."
636 call messages_fatal(1, namespace=namespace)
637 end if
638 if (bitand(how(what_i), option__outputformat__integrate_xy) /= 0) then
639 message(1) = "OutputFormat = integrate_xy not available with Dimensions <= 2."
640 call messages_fatal(1, namespace=namespace)
641 end if
642 if (bitand(how(what_i), option__outputformat__integrate_xz) /= 0) then
643 message(1) = "OutputFormat = integrate_xz not available with Dimensions <= 2."
644 call messages_fatal(1, namespace=namespace)
645 end if
646 if (bitand(how(what_i), option__outputformat__integrate_yz) /= 0) then
647 message(1) = "OutputFormat = integrate_yz not available with Dimensions <= 2."
648 call messages_fatal(1, namespace=namespace)
649 end if
650 if (bitand(how(what_i), option__outputformat__dx) /= 0) then
651 message(1) = "OutputFormat = dx not available with Dimensions <= 2."
652 call messages_fatal(1, namespace=namespace)
653 end if
654 if (bitand(how(what_i), option__outputformat__cube) /= 0) then
655 message(1) = "OutputFormat = cube not available with Dimensions <= 2."
656 call messages_fatal(1, namespace=namespace)
657 end if
658 end if
659
660#if !defined(HAVE_NETCDF)
661 if (bitand(how(what_i), option__outputformat__netcdf) /= 0) then
662 message(1) = 'Octopus was compiled without NetCDF support.'
663 message(2) = 'It is not possible to write output in NetCDF format.'
664 call messages_fatal(2, namespace=namespace)
665 end if
666#endif
667#if !defined(HAVE_ETSF_IO)
668 if (bitand(how(what_i), option__outputformat__etsf) /= 0) then
669 message(1) = 'Octopus was compiled without ETSF_IO support.'
670 message(2) = 'It is not possible to write output in ETSF format.'
671 call messages_fatal(2, namespace=namespace)
672 end if
673#endif
674
675
676 end if
677 end if
678
679 if (output_interval(what_i) < 0) then
680 message(1) = "OutputInterval must be >= 0."
681 call messages_fatal(1, namespace=namespace)
682 end if
683 end do
684
686 end subroutine io_function_read_what_how_when
687
688 ! -------------------------------------------------------------------
691 ! ".", "func", mesh, sb, func, M_ONE, ierr)
692 ! -------------------------------------------------------------------
693 integer(int64) function io_function_fill_how(where) result(how)
694 character(len=*), intent(in) :: where
695
696 push_sub(io_function_fill_how)
697
698 how = 0
699 if (index(where, "AxisX") /= 0) how = ior(how, option__outputformat__axis_x)
700 if (index(where, "AxisY") /= 0) how = ior(how, option__outputformat__axis_y)
701 if (index(where, "AxisZ") /= 0) how = ior(how, option__outputformat__axis_z)
702 if (index(where, "PlaneX") /= 0) how = ior(how, option__outputformat__plane_x)
703 if (index(where, "PlaneY") /= 0) how = ior(how, option__outputformat__plane_y)
704 if (index(where, "PlaneZ") /= 0) how = ior(how, option__outputformat__plane_z)
705 if (index(where, "IntegrateXY") /= 0) how = ior(how, option__outputformat__integrate_xy)
706 if (index(where, "IntegrateXZ") /= 0) how = ior(how, option__outputformat__integrate_xz)
707 if (index(where, "IntegrateYZ") /= 0) how = ior(how, option__outputformat__integrate_yz)
708 if (index(where, "PlaneZ") /= 0) how = ior(how, option__outputformat__plane_z)
709 if (index(where, "DX") /= 0) how = ior(how, option__outputformat__dx)
710 if (index(where, "XCrySDen") /= 0) how = ior(how, option__outputformat__xcrysden)
711 if (index(where, "Binary") /= 0) how = ior(how, option__outputformat__binary)
712 if (index(where, "MeshIndex") /= 0) how = ior(how, option__outputformat__mesh_index)
713 if (index(where, "XYZ") /= 0) how = ior(how, option__outputformat__xyz)
714#if defined(HAVE_NETCDF)
715 if (index(where, "NETCDF") /= 0) how = ior(how, option__outputformat__netcdf)
716#endif
717 if (index(where, "Cube") /= 0) how = ior(how, option__outputformat__cube)
718 if (index(where, "VTK") /= 0) how = ior(how, option__outputformat__vtk)
719
720 pop_sub(io_function_fill_how)
721 end function io_function_fill_how
722
723 ! ---------------------------------------------------------
728 subroutine write_canonicalized_xyz_file(dir, fname, space, latt, pos, atoms, box, namespace)
729 character(len=*), intent(in) :: dir
730 character(len=*), intent(in) :: fname
731 class(space_t), intent(in) :: space
732 type(lattice_vectors_t), intent(in) :: latt
733 real(real64), intent(in) :: pos(:,:)
734 type(atom_t), intent(in) :: atoms(:)
735 class(box_t), intent(in) :: box
736 type(namespace_t), intent(in) :: namespace
737
738 integer :: iunit
739 integer :: idir
740 real(real64) :: position
741 integer :: iatom
742
744
745 call io_mkdir(dir, namespace)
746 iunit = io_open(trim(dir)//'/'//trim(fname)//'.xyz', namespace, action='write', position='asis')
747
748 write(iunit, '(i6)') size(atoms)
749 write(iunit, '(a,a,a)', advance='no') trim(space%short_info()), '; ', trim(box%short_info(unit_angstrom))
750 if (space%is_periodic()) then
751 write(iunit, '(a,a)') '; ', trim(latt%short_info(unit_angstrom))
752 else
753 write(iunit, '()')
754 end if
755
756 ! xyz-style labels and positions:
757 do iatom = 1, size(atoms)
758 write(iunit, '(10a)', advance='no') atoms(iatom)%label
759 do idir = 1, 3
760 if (idir <= space%dim) then
761 position = pos(idir, iatom)
762 else
763 position = m_zero
764 end if
765 write(iunit, '(1x, f11.6)', advance='no') units_from_atomic(unit_angstrom, position)
766 end do
767 write(iunit, '()')
768 end do
769
770 call io_close(iunit)
771
773 end subroutine write_canonicalized_xyz_file
774
775 ! ---------------------------------------------------------
776 subroutine write_xsf_geometry_file(dir, fname, space, latt, pos, atoms, mesh, namespace, total_forces)
777 character(len=*), intent(in) :: dir, fname
778 class(space_t), intent(in) :: space
779 type(lattice_vectors_t), intent(in) :: latt
780 real(real64), intent(in) :: pos(:,:)
781 type(atom_t), intent(in) :: atoms(:)
782 class(mesh_t), intent(in) :: mesh
783 type(namespace_t), intent(in) :: namespace
784 real(real64), optional, intent(in) :: total_forces(:,:)
785
786 integer :: iunit
787 real(real64), allocatable :: forces(:,:)
788
789 if (.not. mpi_grp_is_root(mpi_world)) return
790
792
793 call io_mkdir(dir, namespace)
794 iunit = io_open(trim(dir)//'/'//trim(fname)//'.xsf', namespace, action='write', position='asis')
795
796 if (present(total_forces)) then
797 safe_allocate(forces(1:space%dim, 1:size(atoms)))
798 forces = units_from_atomic(units_out%force, total_forces)
799 call write_xsf_geometry(iunit, space, latt, pos, atoms, mesh, forces = forces)
800 safe_deallocate_a(forces)
801 else
802 call write_xsf_geometry(iunit, space, latt, pos, atoms, mesh)
803 end if
804
805 call io_close(iunit)
806
808 end subroutine write_xsf_geometry_file
809
810 ! ---------------------------------------------------------
813 subroutine write_xsf_geometry(iunit, space, latt, pos, atoms, mesh, forces, index)
814 integer, intent(in) :: iunit
815 class(space_t), intent(in) :: space
816 type(lattice_vectors_t), intent(in) :: latt
817 real(real64), intent(in) :: pos(:,:)
818 type(atom_t), intent(in) :: atoms(:)
819 class(mesh_t), intent(in) :: mesh
820 real(real64), optional, intent(in) :: forces(:, :)
821 integer, optional, intent(in) :: index
822
823 integer :: idir, idir2, iatom, index_, natoms
824 character(len=7) :: index_str
825 real(real64) :: offset(space%dim)
826 real(real64) :: rlattice(space%dim, space%dim)
827
828 push_sub(write_xsf_geometry)
829
830 if (present(index)) then
831 write(index_str, '(a,i6)') ' ', index
832 index_ = index
833 else
834 write(index_str, '(a)') ''
835 index_ = 1
836 end if
837 natoms = size(pos, dim=2)
838
839 ! The corner of the cell is always (0,0,0) to XCrySDen
840 ! so the offset is applied to the atomic coordinates.
841 ! Along periodic dimensions the offset is -1/2 in reduced coordinates, as
842 ! our origin is at the center of the cell instead of being at the edge.
843 offset(1:space%dim) = latt%red_to_cart(spread(-m_half, 1, space%dim))
844 ! Offset in aperiodic directions:
845 do idir = space%periodic_dim + 1, space%dim
846 offset(idir) = -(mesh%idx%ll(idir) - 1)/2 * mesh%spacing(idir)
847 end do
848
849 if (space%is_periodic()) then
850 if (index_ == 1) then
851 select case (space%periodic_dim)
852 case (3)
853 write(iunit, '(a)') 'CRYSTAL'
854 case (2)
855 write(iunit, '(a)') 'SLAB'
856 case (1)
857 write(iunit, '(a)') 'POLYMER'
858 end select
859 end if
860
861 write(iunit, '(a)') 'PRIMVEC'//trim(index_str)
862
863 !Computes the rlattice corresponding to the 3D periodic version of the simulation box
864 rlattice = latt%rlattice
865 do idir = space%periodic_dim+1, space%dim
866 rlattice(:,idir) = rlattice(:,idir)*m_two*mesh%box%bounding_box_l(idir)
867 end do
868
869 do idir = 1, space%dim
870 write(iunit, '(3f12.6)') (units_from_atomic(units_out%length, rlattice(idir2, idir)), idir2 = 1, space%dim)
871 end do
872
873 write(iunit, '(a)') 'PRIMCOORD'//trim(index_str)
874 write(iunit, '(i10, a)') natoms, ' 1'
875 else
876 write(iunit, '(a)') 'ATOMS'//trim(index_str)
877 end if
878
879 ! BoxOffset should be considered here
880 do iatom = 1, natoms
881 write(iunit, '(a10, 3f12.6)', advance='no') trim(atoms(iatom)%label), &
882 (units_from_atomic(units_out%length, pos(idir, iatom) - offset(idir)), idir = 1, space%dim)
883 if (present(forces)) then
884 write(iunit, '(5x, 3f12.6)', advance='no') forces(:, iatom)
885 end if
886 write(iunit, '()')
887 end do
888
889 pop_sub(write_xsf_geometry)
890 end subroutine write_xsf_geometry
891
892! ---------------------------------------------------------
895 subroutine write_xsf_geometry_supercell(iunit, space, latt, pos, atoms, mesh, centers, supercell, extra_atom)
896 integer, intent(in) :: iunit
897 class(space_t), intent(in) :: space
898 type(lattice_vectors_t), intent(in) :: latt
899 real(real64), intent(in) :: pos(:,:)
900 type(atom_t), intent(in) :: atoms(:)
901 class(mesh_t), intent(in) :: mesh
902 real(real64), intent(in) :: centers(:, :)
903 integer, intent(in) :: supercell(:)
904 real(real64), optional, intent(in) :: extra_atom(:)
905
906 integer :: idir, idir2, iatom, index_
907 character(len=7) :: index_str
908 real(real64) :: offset(3)
909 integer :: irep, nreplica
910
912
913 write(index_str, '(a)') ''
914 index_ = 1
915
916 nreplica = product(supercell(1:space%dim))
917
918 ! The corner of the cell is always (0,0,0) to XCrySDen
919 ! so the offset is applied to the atomic coordinates.
920 ! Offset in periodic directions:
921 offset(1:space%dim) = latt%red_to_cart(spread(-m_half, 1, space%dim))
922 offset(1:3) = offset(1:3) + centers(1:3,1)
923 ! Offset in aperiodic directions:
924 do idir = space%periodic_dim + 1, 3
925 offset(idir) = -(mesh%idx%ll(idir) - 1)/2 * mesh%spacing(idir)
926 end do
927
928 if(space%is_periodic()) then
929 if(index_ == 1) then
930 select case(space%periodic_dim)
931 case(3)
932 write(iunit, '(a)') 'CRYSTAL'
933 case(2)
934 write(iunit, '(a)') 'SLAB'
935 case(1)
936 write(iunit, '(a)') 'POLYMER'
937 end select
938 end if
939
940 write(iunit, '(a)') 'PRIMVEC'//trim(index_str)
941
942 do idir = 1, space%dim
943 write(iunit, '(3f12.6)') (units_from_atomic(units_out%length, &
944 latt%rlattice(idir2, idir)*supercell(idir)), idir2 = 1, space%dim)
945 end do
946
947 write(iunit, '(a)') 'PRIMCOORD'//trim(index_str)
948 if(.not.present(extra_atom)) then
949 write(iunit, '(i10, a)') size(atoms)*nreplica, ' 1'
950 else
951 write(iunit, '(i10, a)') size(atoms)*nreplica+1, ' 1'
952 end if
953 else
954 write(iunit, '(a)') 'ATOMS'//trim(index_str)
955 end if
956
957
958 do irep = 1, nreplica
959 ! BoxOffset should be considered here
960 do iatom = 1, size(atoms)
961 write(iunit, '(a10, 3f12.6)', advance='no') trim(atoms(iatom)%label), &
962 (units_from_atomic(units_out%length, pos(idir, iatom) + centers(idir, irep) &
963 - offset(idir)), idir = 1, space%dim)
964 write(iunit, '()')
965 end do
966 end do
967 write(iunit, '(a10, 3f12.6)', advance='no') 'X', &
968 (units_from_atomic(units_out%length, extra_atom(idir) - offset(idir)), idir = 1, space%dim)
969 write(iunit, '()')
970
973
974
975#if defined(HAVE_NETCDF)
976 ! ---------------------------------------------------------
977 subroutine ncdf_error(func, status, filename, namespace, ierr)
978 character(len=*), intent(in) :: func
979 integer, intent(in) :: status
980 character(len=*), intent(in) :: filename
981 type(namespace_t), intent(in) :: namespace
982 integer, intent(inout) :: ierr
983
984 push_sub(ncdf_error)
985
986 if (status == nf90_noerr) then
987 pop_sub(ncdf_error)
988 return
989 end if
990
991 write(message(1),'(3a)') "NETCDF error in function '" , trim(func) , "'"
992 write(message(2),'(3a)') "(reading/writing ", trim(filename) , ")"
993 write(message(3), '(6x,a,a)')'Error code = ', trim(nf90_strerror(status))
994 call messages_warning(3, namespace=namespace)
995 ierr = 5
996
997 pop_sub(ncdf_error)
998 end subroutine ncdf_error
999#endif
1000
1001 ! ---------------------------------------------------------
1002 subroutine transpose3(in, out)
1003 real(real64), intent(in) :: in(:, :, :)
1004 real(real64), intent(out) :: out(:, :, :)
1005 integer :: ix, iy, iz
1006
1007 push_sub(transpose3)
1008
1009 do ix = lbound(in, 1), ubound(in, 1)
1010 do iy = lbound(in, 2), ubound(in, 2)
1011 do iz = lbound(in, 3), ubound(in, 3)
1012 out(iz, iy, ix) = in(ix, iy, iz)
1013 end do
1014 end do
1015 end do
1016
1017 pop_sub(transpose3)
1018 end subroutine transpose3
1019
1020#include "undef.F90"
1021#include "real.F90"
1022#include "io_function_inc.F90"
1023
1024#include "undef.F90"
1025#include "complex.F90"
1026#include "io_function_inc.F90"
1027#include "undef.F90"
1028
1029end module io_function_oct_m
1030
1031!! Local Variables:
1032!! mode: f90
1033!! coding: utf-8
1034!! End:
Fast Fourier Transform module. This module provides a single interface that works with different FFT ...
Definition: fft.F90:118
real(real64), parameter, public m_two
Definition: global.F90:190
real(real64), parameter, public m_zero
Definition: global.F90:188
real(real64), parameter, public m_half
Definition: global.F90:194
This module implements the index, used for the mesh points.
Definition: index.F90:122
subroutine zio_function_output_vector_bz(how, dir, fname, namespace, space, kpt, kpoints, ff, unit, ierr, grp, root)
subroutine, public write_canonicalized_xyz_file(dir, fname, space, latt, pos, atoms, box, namespace)
Write canonicalized xyz file with atom labels and positions in Angstroms. Includes information about ...
subroutine, public zio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
subroutine, public dio_function_output_global(how, dir, fname, namespace, space, mesh, ff, unit, ierr)
subroutine, public io_function_read_what_how_when(namespace, space, what, how, output_interval, what_tag_in, how_tag_in, output_interval_tag_in, ignore_error)
subroutine, public zio_function_input(filename, namespace, space, mesh, ff, ierr, map)
Reads a mesh function from file filename, and puts it into ff. If the map argument is passed,...
subroutine, public write_xsf_geometry(iunit, space, latt, pos, atoms, mesh, forces, index)
for format specification see: http:
subroutine zio_function_output_vector(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
subroutine write_xsf_geometry_supercell(iunit, space, latt, pos, atoms, mesh, centers, supercell, extra_atom)
for format specification see: http:
subroutine dio_function_output_supercell(how, dir, fname, mesh, space, latt, ff, centers, supercell, unit, ierr, namespace, pos, atoms, grp, root, is_global, extra_atom)
subroutine, public dio_function_output(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
Top-level IO routine for functions defined on the mesh.
subroutine dio_function_output_vector(how, dir, fname, namespace, space, mesh, ff, unit, ierr, pos, atoms, grp, root)
subroutine transpose3(in, out)
subroutine zio_function_output_supercell(how, dir, fname, mesh, space, latt, ff, centers, supercell, unit, ierr, namespace, pos, atoms, grp, root, is_global, extra_atom)
subroutine, public write_xsf_geometry_file(dir, fname, space, latt, pos, atoms, mesh, namespace, total_forces)
subroutine, public dout_cf_netcdf(filename, ierr, cf, cube, space, spacing, transpose, unit, namespace)
Writes a cube_function in netcdf format.
subroutine dio_function_output_global_bz(how, dir, fname, namespace, kpoints, ff, unit, ierr)
subroutine ncdf_error(func, status, filename, namespace, ierr)
integer, parameter, private zoutput_kind
subroutine, public zio_function_output_global(how, dir, fname, namespace, space, mesh, ff, unit, ierr)
subroutine, public zout_cf_netcdf(filename, ierr, cf, cube, space, spacing, transpose, unit, namespace)
Writes a cube_function in netcdf format.
subroutine dio_function_output_vector_bz(how, dir, fname, namespace, space, kpt, kpoints, ff, unit, ierr, grp, root)
integer(int64) function, public io_function_fill_how(where)
Use this function to quickly plot functions for debugging purposes: call dio_function_output(io_funct...
subroutine, public dio_function_input(filename, namespace, space, mesh, ff, ierr, map)
Reads a mesh function from file filename, and puts it into ff. If the map argument is passed,...
subroutine zio_function_output_global_bz(how, dir, fname, namespace, kpoints, ff, unit, ierr)
Definition: io.F90:114
subroutine, public io_close(iunit, grp)
Definition: io.F90:418
subroutine, public io_mkdir(fname, namespace, parents)
Definition: io.F90:311
integer function, public io_open(file, namespace, action, status, form, position, die, recl, grp)
Definition: io.F90:352
This module defines various routines, operating on mesh functions.
This module defines the meshes, which are used in Octopus.
Definition: mesh.F90:118
subroutine, public messages_variable_is_block(namespace, name)
Definition: messages.F90:1069
subroutine, public messages_warning(no_lines, all_nodes, namespace)
Definition: messages.F90:537
subroutine, public messages_obsolete_variable(namespace, name, rep)
Definition: messages.F90:1045
character(len=256), dimension(max_lines), public message
to be output by fatal, warning
Definition: messages.F90:160
subroutine, public messages_fatal(no_lines, only_root_writes, namespace)
Definition: messages.F90:414
subroutine, public messages_input_error(namespace, var, details, row, column)
Definition: messages.F90:713
logical function mpi_grp_is_root(grp)
Is the current MPI process of grpcomm, root.
Definition: mpi.F90:430
type(mpi_grp_t), public mpi_world
Definition: mpi.F90:266
Some general things and nomenclature:
Definition: par_vec.F90:171
integer function, public parse_block(namespace, name, blk, check_varinfo_)
Definition: parser.F90:618
brief This module defines the class unit_t which is used by the unit_systems_oct_m module.
Definition: unit.F90:132
This module defines the unit system, used for input and output.
type(unit_t), public unit_angstrom
For XYZ files.
type(unit_system_t), public units_out
This module is intended to contain simple general-purpose utility functions and procedures.
Definition: utils.F90:118
class to tell whether a point is inside or outside
Definition: box.F90:141
Describes mesh distribution to nodes.
Definition: mesh.F90:186
int true(void)