Properly finalize MPI on mdrun -version. Fixes #1313
[gromacs.git] / include / sfactor.h
blob7f6d9bc4a0a4dc0b58b6a455287e7514eded1b1a
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 #ifndef _sfactor_h
40 #define _sfactor_h
42 #include "visibility.h"
43 #include "index.h"
44 #include "types/simple.h"
45 #include "gmxcomplex.h"
46 #include "oenv.h"
50 #ifdef __cplusplus
51 extern "C" {
52 #endif
55 typedef struct gmx_structurefactors gmx_structurefactors_t;
57 typedef struct structure_factor structure_factor_t;
59 typedef struct reduced_atom reduced_atom_t;
61 int * create_indexed_atom_type (reduced_atom_t * atm, int size);
63 void compute_structure_factor (structure_factor_t * sft, matrix box,
64 reduced_atom_t * red, int isize, real start_q,
65 real end_q, int group, real **sf_table);
67 gmx_structurefactors_t *gmx_structurefactors_init(const char *datfn);
69 void gmx_structurefactors_done(gmx_structurefactors_t *gsf);
71 int gmx_structurefactors_get_sf(gmx_structurefactors_t *gsf, int elem, real a[4], real b[4], real *c);
73 real **gmx_structurefactors_table(gmx_structurefactors_t *gsf, real momentum, real ref_k,
74 real lambda, int n_angles);
76 void save_data (structure_factor_t * sft, const char *file, int ngrps,
77 real start_q, real end_q, const output_env_t oenv);
79 double CMSF (gmx_structurefactors_t *gsf, int type, int nh, double lambda, double sin_theta);
81 int return_atom_type (const char *name, gmx_structurefactors_t *gsf);
83 void rearrange_atoms (reduced_atom_t * positions, t_trxframe *fr, atom_id * index,
84 int isize, t_topology * top, gmx_bool flag, gmx_structurefactors_t *gsf);
86 GMX_LIBGMX_EXPORT
87 int do_scattering_intensity (const char* fnTPS, const char* fnNDX,
88 const char* fnXVG, const char *fnTRX,
89 const char* fnDAT,
90 real start_q, real end_q,
91 real energy, int ng, const output_env_t oenv);
93 t_complex *** rc_tensor_allocation(int x, int y, int z);
95 real **compute_scattering_factor_table (gmx_structurefactors_t *gsf, structure_factor_t * sft, int *nsftable);
97 #ifdef __cplusplus
99 #endif
100 #endif