Work around memory leak in t_state
[gromacs.git] / src / gromacs / fileio / xvgr.cpp
blob6034a4f811eaa2e3f3cc8155647d2f9e9d8a494d
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 * Copyright (c) 2013,2014,2015, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "xvgr.h"
41 #include <cctype>
42 #include <cstring>
44 #include <string>
46 #include "gromacs/fileio/gmxfio.h"
47 #include "gromacs/fileio/oenv.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/utility/binaryinformation.h"
50 #include "gromacs/utility/coolstuff.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/exceptions.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
56 #include "gromacs/utility/sysinfo.h"
58 gmx_bool output_env_get_print_xvgr_codes(const gmx_output_env_t *oenv)
60 int xvg_format;
62 xvg_format = output_env_get_xvg_format(oenv);
64 return (xvg_format == exvgXMGRACE || xvg_format == exvgXMGR);
67 static char *xvgrstr(const char *gmx, const gmx_output_env_t *oenv,
68 char *buf, int buflen)
70 /* Supported greek letter names and corresponding xmgrace/xmgr symbols */
71 const char *sym[] = { "beta", "chi", "delta", "eta", "lambda", "mu", "omega", "phi", "psi", "rho", "theta", NULL };
72 const char symc[] = { 'b', 'c', 'd', 'h', 'l', 'm', 'w', 'f', 'y', 'r', 'q', '\0' };
73 int xvgf;
74 gmx_bool bXVGR;
75 int g, b, i;
76 char c;
78 xvgf = output_env_get_xvg_format(oenv);
79 bXVGR = (xvgf == exvgXMGRACE || xvgf == exvgXMGR);
81 g = 0;
82 b = 0;
83 while (gmx[g] != '\0')
85 /* Check with the largest string we have ("lambda"), add one for \0 */
86 if (b + 6 + 1 >= buflen)
88 gmx_fatal(FARGS, "Output buffer length in xvgstr (%d) too small to process xvg input string '%s'", buflen, gmx);
90 if (gmx[g] == '\\')
92 g++;
93 if (gmx[g] == 's')
95 /* Subscript */
96 if (bXVGR)
98 buf[b++] = '\\';
99 buf[b++] = 's';
101 else
103 buf[b++] = '_';
105 g++;
107 else if (gmx[g] == 'S')
109 /* Superscript */
110 if (bXVGR)
112 buf[b++] = '\\';
113 buf[b++] = 'S';
115 else
117 buf[b++] = '^';
119 g++;
121 else if (gmx[g] == 'N')
123 /* End sub/superscript */
124 if (bXVGR)
126 buf[b++] = '\\';
127 buf[b++] = 'N';
129 else
131 if (gmx[g+1] != ' ')
133 buf[b++] = ' ';
136 g++;
138 else if (gmx[g] == '4')
140 /* Backward compatibility for xmgr normal font "\4" */
141 switch (xvgf)
143 case exvgXMGRACE:
144 sprintf(buf+b, "%s", "\\f{}");
145 break;
146 case exvgXMGR:
147 sprintf(buf+b, "%s", "\\4");
148 break;
149 default:
150 buf[b] = '\0';
151 break;
153 g++;
154 b = strlen(buf);
156 else if (gmx[g] == '8')
158 /* Backward compatibility for xmgr symbol font "\8" */
159 switch (xvgf)
161 case exvgXMGRACE:
162 sprintf(buf+b, "%s", "\\x");
163 break;
164 case exvgXMGR:
165 sprintf(buf+b, "%s", "\\8");
166 break;
167 default:
168 buf[b] = '\0';
169 break;
171 g++;
172 b = std::strlen(buf);
174 else
176 /* Check for special symbol */
177 i = 0;
178 while (sym[i] != NULL &&
179 gmx_strncasecmp(sym[i], gmx+g, std::strlen(sym[i])) != 0)
181 i++;
183 if (sym[i] != NULL)
185 c = symc[i];
186 if (std::isupper(gmx[g]))
188 c = std::toupper(c);
190 switch (xvgf)
192 case exvgXMGRACE:
193 sprintf(buf+b, "%s%c%s", "\\x", c, "\\f{}");
194 break;
195 case exvgXMGR:
196 sprintf(buf+b, "%s%c%s", "\\8", c, "\\4");
197 break;
198 default:
199 std::strncat(buf+b, gmx+g, std::strlen(sym[i]));
200 b += std::strlen(sym[i]);
201 if (gmx[g+std::strlen(sym[i])] != ' ')
203 buf[b++] = ' ';
205 buf[b] = '\0';
206 break;
208 g += std::strlen(sym[i]);
209 b = std::strlen(buf);
211 else
213 /* Unknown escape sequence, this should not happen.
214 * We do not generate a fatal error, since that might
215 * stop otherwise functioning code from working.
216 * Copy the backslash to the output and continue processing.
218 buf[b++] = '\\';
222 else
224 buf[b++] = gmx[g++];
228 buf[b++] = '\0';
230 return buf;
233 void xvgr_header(FILE *fp, const char *title, const char *xaxis,
234 const char *yaxis, int exvg_graph_type,
235 const gmx_output_env_t *oenv)
237 char buf[STRLEN];
239 if (output_env_get_print_xvgr_codes(oenv))
241 gmx_format_current_time(buf, STRLEN);
242 fprintf(fp, "# This file was created %s", buf);
245 gmx::BinaryInformationSettings settings;
246 settings.generatedByHeader(true);
247 settings.linePrefix("# ");
248 gmx::printBinaryInformation(fp, output_env_get_program_context(oenv),
249 settings);
251 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
252 fprintf(fp, "# %s is part of G R O M A C S:\n#\n",
253 output_env_get_program_display_name(oenv));
254 fprintf(fp, "# %s\n#\n", gmx::bromacs().c_str());
255 fprintf(fp, "@ title \"%s\"\n", xvgrstr(title, oenv, buf, STRLEN));
256 fprintf(fp, "@ xaxis label \"%s\"\n",
257 xvgrstr(xaxis, oenv, buf, STRLEN));
258 fprintf(fp, "@ yaxis label \"%s\"\n",
259 xvgrstr(yaxis, oenv, buf, STRLEN));
260 switch (exvg_graph_type)
262 case exvggtXNY:
263 if (output_env_get_xvg_format(oenv) == exvgXMGR)
265 fprintf(fp, "@TYPE nxy\n");
267 else
269 fprintf(fp, "@TYPE xy\n");
271 break;
272 case exvggtXYDY:
273 fprintf(fp, "@TYPE xydy\n");
274 break;
275 case exvggtXYDYDY:
276 fprintf(fp, "@TYPE xydydy\n");
277 break;
282 FILE *xvgropen_type(const char *fn, const char *title, const char *xaxis,
283 const char *yaxis, int exvg_graph_type,
284 const gmx_output_env_t *oenv)
286 FILE *fp;
288 fp = gmx_fio_fopen(fn, "w");
290 xvgr_header(fp, title, xaxis, yaxis, exvg_graph_type, oenv);
292 return fp;
295 FILE *xvgropen(const char *fn, const char *title, const char *xaxis,
296 const char *yaxis, const gmx_output_env_t *oenv)
298 return xvgropen_type(fn, title, xaxis, yaxis, exvggtXNY, oenv);
301 void
302 xvgrclose(FILE *fp)
304 gmx_fio_fclose(fp);
307 void xvgr_subtitle(FILE *out, const char *subtitle, const gmx_output_env_t *oenv)
309 char buf[STRLEN];
311 if (output_env_get_print_xvgr_codes(oenv))
313 fprintf(out, "@ subtitle \"%s\"\n", xvgrstr(subtitle, oenv, buf, STRLEN));
317 void xvgr_view(FILE *out, real xmin, real ymin, real xmax, real ymax,
318 const gmx_output_env_t *oenv)
320 if (output_env_get_print_xvgr_codes(oenv))
322 fprintf(out, "@ view %g, %g, %g, %g\n", xmin, ymin, xmax, ymax);
326 void xvgr_world(FILE *out, real xmin, real ymin, real xmax, real ymax,
327 const gmx_output_env_t *oenv)
329 if (output_env_get_print_xvgr_codes(oenv))
331 fprintf(out, "@ world xmin %g\n"
332 "@ world ymin %g\n"
333 "@ world xmax %g\n"
334 "@ world ymax %g\n", xmin, ymin, xmax, ymax);
338 void xvgr_legend(FILE *out, int nsets, const char** setname,
339 const gmx_output_env_t *oenv)
341 int i;
342 char buf[STRLEN];
344 if (output_env_get_print_xvgr_codes(oenv))
346 xvgr_view(out, 0.15, 0.15, 0.75, 0.85, oenv);
347 fprintf(out, "@ legend on\n");
348 fprintf(out, "@ legend box on\n");
349 fprintf(out, "@ legend loctype view\n");
350 fprintf(out, "@ legend %g, %g\n", 0.78, 0.8);
351 fprintf(out, "@ legend length %d\n", 2);
352 for (i = 0; (i < nsets); i++)
354 if (setname[i])
356 if (output_env_get_xvg_format(oenv) == exvgXMGR)
358 fprintf(out, "@ legend string %d \"%s\"\n",
359 i, xvgrstr(setname[i], oenv, buf, STRLEN));
361 else
363 fprintf(out, "@ s%d legend \"%s\"\n",
364 i, xvgrstr(setname[i], oenv, buf, STRLEN));
371 void xvgr_new_dataset(FILE *out, int nr_first, int nsets,
372 const char **setname,
373 const gmx_output_env_t *oenv)
375 int i;
376 char buf[STRLEN];
378 if (output_env_get_print_xvgr_codes(oenv))
380 fprintf(out, "@\n");
381 for (i = 0; (i < nsets); i++)
383 if (setname[i])
385 if (output_env_get_xvg_format(oenv) == exvgXMGR)
387 fprintf(out, "@ legend string %d \"%s\"\n",
388 i+nr_first, xvgrstr(setname[i], oenv, buf, STRLEN));
390 else
392 fprintf(out, "@ s%d legend \"%s\"\n",
393 i+nr_first, xvgrstr(setname[i], oenv, buf, STRLEN));
398 else
400 fprintf(out, "\n");
404 void xvgr_line_props(FILE *out, int NrSet, int LineStyle, int LineColor,
405 const gmx_output_env_t *oenv)
407 if (output_env_get_print_xvgr_codes(oenv))
409 fprintf(out, "@ with g0\n");
410 fprintf(out, "@ s%d linestyle %d\n", NrSet, LineStyle);
411 fprintf(out, "@ s%d color %d\n", NrSet, LineColor);
415 static const char *LocTypeStr[] = { "view", "world" };
416 static const char *BoxFillStr[] = { "none", "color", "pattern" };
418 void xvgr_box(FILE *out,
419 int LocType,
420 real xmin, real ymin, real xmax, real ymax,
421 int LineStyle, int LineWidth, int LineColor,
422 int BoxFill, int BoxColor, int BoxPattern, const gmx_output_env_t *oenv)
424 if (output_env_get_print_xvgr_codes(oenv))
426 fprintf(out, "@with box\n");
427 fprintf(out, "@ box on\n");
428 fprintf(out, "@ box loctype %s\n", LocTypeStr[LocType]);
429 fprintf(out, "@ box %g, %g, %g, %g\n", xmin, ymin, xmax, ymax);
430 fprintf(out, "@ box linestyle %d\n", LineStyle);
431 fprintf(out, "@ box linewidth %d\n", LineWidth);
432 fprintf(out, "@ box color %d\n", LineColor);
433 fprintf(out, "@ box fill %s\n", BoxFillStr[BoxFill]);
434 fprintf(out, "@ box fill color %d\n", BoxColor);
435 fprintf(out, "@ box fill pattern %d\n", BoxPattern);
436 fprintf(out, "@box def\n");
440 /* reads a line into ptr, adjusting len and renewing ptr if neccesary */
441 static char *fgets3(FILE *fp, char **ptr, int *len, int maxlen)
443 int len_remaining = *len; /* remaining amount of allocated bytes in buf */
444 int curp = 0; /* current position in buf to read into */
448 if (len_remaining < 2)
450 if (*len + STRLEN < maxlen)
452 /* This line is longer than len characters, let's increase len! */
453 *len += STRLEN;
454 len_remaining += STRLEN;
455 srenew(*ptr, *len);
457 else
459 /*something is wrong, we'll just keep reading and return NULL*/
460 len_remaining = STRLEN;
461 curp = 0;
464 if (fgets(*ptr + curp, len_remaining, fp) == NULL)
466 /* if last line, skip */
467 return NULL;
469 curp += len_remaining-1; /* overwrite the nul char in next iteration */
470 len_remaining = 1;
472 while ((std::strchr(*ptr, '\n') == NULL) && (!feof(fp)));
474 if (*len + STRLEN >= maxlen)
476 return NULL; /* this line was too long */
479 if (feof(fp))
481 /* We reached EOF before '\n', skip this last line. */
482 return NULL;
485 /* now remove newline */
486 int slen = std::strlen(*ptr);
487 if ((*ptr)[slen-1] == '\n')
489 (*ptr)[slen-1] = '\0';
493 return *ptr;
496 static int wordcount(char *ptr)
498 int i, n = 0, is[2];
499 int cur = 0;
500 #define prev (1-cur)
502 if (NULL != ptr)
504 for (i = 0; (ptr[i] != '\0'); i++)
506 is[cur] = std::isspace(ptr[i]);
507 if ((0 == i) && !is[cur])
509 n++;
511 else if ((i > 0) && (!is[cur] && is[prev]))
513 n++;
515 cur = prev;
518 return n;
521 static char *read_xvgr_string(const char *line)
523 const char *ptr0, *ptr1;
524 char *str;
526 ptr0 = std::strchr(line, '"');
527 if (ptr0 != NULL)
529 ptr0++;
530 ptr1 = std::strchr(ptr0, '"');
531 if (ptr1 != NULL)
533 str = gmx_strdup(ptr0);
534 str[ptr1-ptr0] = '\0';
536 else
538 str = gmx_strdup("");
541 else
543 str = gmx_strdup("");
546 return str;
549 int read_xvg_legend(const char *fn, double ***y, int *ny,
550 char **subtitle, char ***legend)
552 FILE *fp;
553 char *ptr;
554 char *base = NULL;
555 char *fmt = NULL;
556 int k, line = 0, nny, nx, maxx, rval, legend_nalloc, set, nchar;
557 double lf;
558 double **yy = NULL;
559 char *tmpbuf;
560 int len = STRLEN;
561 *ny = 0;
562 nny = 0;
563 nx = 0;
564 maxx = 0;
565 fp = gmx_fio_fopen(fn, "r");
567 snew(tmpbuf, len);
568 if (subtitle != NULL)
570 *subtitle = NULL;
572 legend_nalloc = 0;
573 if (legend != NULL)
575 *legend = NULL;
578 while ((ptr = fgets3(fp, &tmpbuf, &len, 10*STRLEN)) != NULL && ptr[0] != '&')
580 line++;
581 trim(ptr);
582 if (ptr[0] == '@')
584 if (legend != NULL)
586 ptr++;
587 trim(ptr);
588 set = -1;
589 if (std::strncmp(ptr, "subtitle", 8) == 0)
591 ptr += 8;
592 if (subtitle != NULL)
594 *subtitle = read_xvgr_string(ptr);
597 else if (std::strncmp(ptr, "legend string", 13) == 0)
599 ptr += 13;
600 sscanf(ptr, "%d%n", &set, &nchar);
601 ptr += nchar;
603 else if (ptr[0] == 's')
605 ptr++;
606 sscanf(ptr, "%d%n", &set, &nchar);
607 ptr += nchar;
608 trim(ptr);
609 if (std::strncmp(ptr, "legend", 6) == 0)
611 ptr += 6;
613 else
615 set = -1;
618 if (set >= 0)
620 if (set >= legend_nalloc)
622 legend_nalloc = set + 1;
623 srenew(*legend, legend_nalloc);
624 (*legend)[set] = read_xvgr_string(ptr);
629 else if (ptr[0] != '#')
631 if (nny == 0)
633 (*ny) = nny = wordcount(ptr);
634 /* fprintf(stderr,"There are %d columns in your file\n",nny);*/
635 if (nny == 0)
637 return 0;
639 snew(yy, nny);
640 snew(fmt, 3*nny+1);
641 snew(base, 3*nny+1);
643 /* Allocate column space */
644 if (nx >= maxx)
646 maxx += 1024;
647 for (k = 0; (k < nny); k++)
649 srenew(yy[k], maxx);
652 /* Initiate format string */
653 fmt[0] = '\0';
654 base[0] = '\0';
656 /* fprintf(stderr,"ptr='%s'\n",ptr);*/
657 for (k = 0; (k < nny); k++)
659 std::strcpy(fmt, base);
660 std::strcat(fmt, "%lf");
661 rval = sscanf(ptr, fmt, &lf);
662 /* fprintf(stderr,"rval = %d\n",rval);*/
663 if ((rval == EOF) || (rval == 0))
665 break;
667 yy[k][nx] = lf;
668 srenew(fmt, 3*(nny+1)+1);
669 srenew(base, 3*nny+1);
670 std::strcat(base, "%*s");
672 if (k != nny)
674 fprintf(stderr, "Only %d columns on line %d in file %s\n",
675 k, line, fn);
676 for (; (k < nny); k++)
678 yy[k][nx] = 0.0;
681 nx++;
684 gmx_fio_fclose(fp);
686 *y = yy;
687 sfree(tmpbuf);
688 sfree(base);
689 sfree(fmt);
691 if (legend_nalloc > 0)
693 if (*ny - 1 > legend_nalloc)
695 srenew(*legend, *ny-1);
696 for (set = legend_nalloc; set < *ny-1; set++)
698 (*legend)[set] = NULL;
703 return nx;
706 int read_xvg(const char *fn, double ***y, int *ny)
708 return read_xvg_legend(fn, y, ny, NULL, NULL);
711 void write_xvg(const char *fn, const char *title, int nx, int ny, real **y,
712 const char **leg, const gmx_output_env_t *oenv)
714 FILE *fp;
715 int i, j;
717 fp = xvgropen(fn, title, "X", "Y", oenv);
718 if (leg)
720 xvgr_legend(fp, ny-1, leg, oenv);
722 for (i = 0; (i < nx); i++)
724 for (j = 0; (j < ny); j++)
726 fprintf(fp, " %12.5e", y[j][i]);
728 fprintf(fp, "\n");
730 xvgrclose(fp);
733 real **read_xvg_time(const char *fn,
734 gmx_bool bHaveT, gmx_bool bTB, real tb, gmx_bool bTE, real te,
735 int nsets_in, int *nset, int *nval, real *dt, real **t)
737 FILE *fp;
738 #define MAXLINELEN 16384
739 char line0[MAXLINELEN];
740 char *line;
741 int t_nalloc, *val_nalloc, a, narg, n, sin, set, nchar;
742 double dbl;
743 gmx_bool bEndOfSet, bTimeInRange, bFirstLine = TRUE;
744 real **val;
746 t_nalloc = 0;
747 *t = NULL;
748 val = NULL;
749 val_nalloc = NULL;
750 *nset = 0;
751 *dt = 0;
752 fp = gmx_fio_fopen(fn, "r");
753 for (sin = 0; sin < nsets_in; sin++)
755 if (nsets_in == 1)
757 narg = 0;
759 else
761 narg = bHaveT ? 2 : 1;
763 n = 0;
764 bEndOfSet = FALSE;
765 while (!bEndOfSet && fgets(line0, MAXLINELEN, fp))
767 line = line0;
768 /* Remove whitespace */
769 while (line[0] == ' ' || line[0] == '\t')
771 line++;
773 bEndOfSet = (line[0] == '&');
774 if (line[0] != '#' && line[0] != '@' && line[0] != '\n' && !bEndOfSet)
776 if (bFirstLine && bHaveT)
778 /* Check the first line that should contain data */
779 a = sscanf(line, "%lf%lf", &dbl, &dbl);
780 if (a == 0)
782 gmx_fatal(FARGS, "Expected a number in %s on line:\n%s", fn, line0);
784 else if (a == 1)
786 fprintf(stderr, "Found only 1 number on line, "
787 "assuming no time is present.\n");
788 bHaveT = FALSE;
789 if (nsets_in > 1)
791 narg = 1;
796 a = 0;
797 bTimeInRange = TRUE;
798 while ((a < narg || (nsets_in == 1 && n == 0)) &&
799 sscanf(line, "%lf%n", &dbl, &nchar) == 1 && bTimeInRange)
801 /* Use set=-1 as the time "set" */
802 if (sin)
804 if (!bHaveT || (a > 0))
806 set = sin;
808 else
810 set = -1;
813 else
815 if (!bHaveT)
817 set = a;
819 else
821 set = a-1;
824 if (set == -1 && ((bTB && dbl < tb) || (bTE && dbl > te)))
826 bTimeInRange = FALSE;
829 if (bTimeInRange)
831 if (n == 0)
833 if (nsets_in == 1)
835 narg++;
837 if (set >= 0)
839 *nset = set+1;
840 srenew(val, *nset);
841 srenew(val_nalloc, *nset);
842 val_nalloc[set] = 0;
843 val[set] = NULL;
846 if (set == -1)
848 if (sin == 0)
850 if (n >= t_nalloc)
852 t_nalloc = over_alloc_small(n);
853 srenew(*t, t_nalloc);
855 (*t)[n] = dbl;
857 /* else we should check the time of the next sets with set 0 */
859 else
861 if (n >= val_nalloc[set])
863 val_nalloc[set] = over_alloc_small(n);
864 srenew(val[set], val_nalloc[set]);
866 val[set][n] = (real)dbl;
869 a++;
870 line += nchar;
872 if (line0[strlen(line0)-1] != '\n')
874 fprintf(stderr, "File %s does not end with a newline, ignoring the last line\n", fn);
876 else if (bTimeInRange)
878 if (a == 0)
880 fprintf(stderr, "Ignoring invalid line in %s:\n%s", fn, line0);
882 else
884 if (a != narg)
886 fprintf(stderr, "Invalid line in %s:\n%s"
887 "Using zeros for the last %d sets\n",
888 fn, line0, narg-a);
890 n++;
893 if (a > 0)
895 bFirstLine = FALSE;
899 if (sin == 0)
901 *nval = n;
902 if (!bHaveT)
904 snew(*t, n);
905 for (a = 0; a < n; a++)
907 (*t)[a] = a;
910 if (n > 1)
912 *dt = (real)((*t)[n-1]-(*t)[0])/(n-1.0);
914 else
916 *dt = 1;
919 else
921 if (n < *nval)
923 fprintf(stderr, "Set %d is shorter (%d) than the previous set (%d)\n",
924 sin+1, n, *nval);
925 *nval = n;
926 fprintf(stderr, "Will use only the first %d points of every set\n",
927 *nval);
931 gmx_fio_fclose(fp);
933 sfree(val_nalloc);
935 return val;