Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
mpim-sw
cdo
Commits
19eb9b95
Commit
19eb9b95
authored
Jul 16, 2011
by
Uwe Schulzweida
Browse files
added Module FC (operators: four2spec spec2four four2grid grid2four)
parent
36970cef
Changes
11
Hide whitespace changes
Inline
Side-by-side
.gitattributes
View file @
19eb9b95
...
...
@@ -291,6 +291,7 @@ src/Eof3d.c -text
src/Eofcoeff.c -text
src/Eofcoeff3d.c -text
src/Exprf.c -text
src/FC.c -text
src/Filedes.c -text
src/Fillmiss.c -text
src/Filter.c -text
...
...
src/FC.c
0 → 100644
View file @
19eb9b95
/*
This file is part of CDO. CDO is a collection of Operators to
manipulate and analyse Climate model Data.
Copyright (C) 2003-2011 Uwe Schulzweida, Uwe.Schulzweida@zmaw.de
See COPYING file for copying and redistribution conditions.
This program is free software; you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation; version 2 of the License.
This program is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
*/
/*
This module contains the following operators:
FC fc2sp Fourier to spectral
FC sp2fc Spectral to fourier
FC fc2gp Fourier to gridpoint
FC gp2fc Gridpoint to fourier
*/
#include
<cdi.h>
#include
"cdo.h"
#include
"cdo_int.h"
#include
"pstream.h"
#include
"specspace.h"
#include
"list.h"
void
*
FC
(
void
*
argument
)
{
int
FC2SP
,
SP2FC
,
FC2GP
,
GP2FC
;
int
operatorID
;
int
streamID1
,
streamID2
;
int
nrecs
,
nvars
;
int
tsID
,
recID
,
varID
,
levelID
;
int
gridsize
;
int
index
,
ngrids
;
int
vlistID1
,
vlistID2
;
int
gridIDsp
=
-
1
,
gridIDgp
=
-
1
,
gridIDfc
=
-
1
;
int
gridID1
=
-
1
,
gridID2
=
-
1
;
int
gridID
;
int
nmiss
;
int
*
vars
;
int
lcopy
=
FALSE
;
int
taxisID1
,
taxisID2
;
int
nlon
=
0
,
nlat
=
0
,
ntr
=
0
;
int
nsp
=
0
,
nfc
=
0
;
double
*
array1
=
NULL
,
*
array2
=
NULL
;
SPTRANS
*
sptrans
=
NULL
;
cdoInitialize
(
argument
);
FC2SP
=
cdoOperatorAdd
(
"fc2sp"
,
0
,
0
,
NULL
);
SP2FC
=
cdoOperatorAdd
(
"sp2fc"
,
0
,
0
,
NULL
);
FC2GP
=
cdoOperatorAdd
(
"fc2gp"
,
0
,
0
,
NULL
);
GP2FC
=
cdoOperatorAdd
(
"gp2fc"
,
0
,
0
,
NULL
);
operatorID
=
cdoOperatorID
();
if
(
UNCHANGED_RECORD
)
lcopy
=
TRUE
;
streamID1
=
streamOpenRead
(
cdoStreamName
(
0
));
vlistID1
=
streamInqVlist
(
streamID1
);
vlistID2
=
vlistDuplicate
(
vlistID1
);
taxisID1
=
vlistInqTaxis
(
vlistID1
);
taxisID2
=
taxisDuplicate
(
taxisID1
);
vlistDefTaxis
(
vlistID2
,
taxisID2
);
ngrids
=
vlistNgrids
(
vlistID1
);
/* find first spectral grid */
for
(
index
=
0
;
index
<
ngrids
;
index
++
)
{
gridID
=
vlistGrid
(
vlistID1
,
index
);
if
(
gridInqType
(
gridID
)
==
GRID_SPECTRAL
)
{
gridIDsp
=
gridID
;
break
;
}
}
/* find first gaussian grid */
for
(
index
=
0
;
index
<
ngrids
;
index
++
)
{
gridID
=
vlistGrid
(
vlistID1
,
index
);
if
(
gridInqType
(
gridID
)
==
GRID_GAUSSIAN
)
{
gridIDgp
=
gridID
;
break
;
}
}
/* find first fourier grid */
for
(
index
=
0
;
index
<
ngrids
;
index
++
)
{
gridID
=
vlistGrid
(
vlistID1
,
index
);
if
(
gridInqType
(
gridID
)
==
GRID_FOURIER
)
{
gridIDfc
=
gridID
;
break
;
}
}
/* define output grid */
if
(
operatorID
==
FC2SP
)
{
if
(
gridIDfc
==
-
1
)
cdoWarning
(
"No fourier data found!"
);
gridID1
=
gridIDfc
;
if
(
gridID1
!=
-
1
)
{
nfc
=
gridInqSize
(
gridID1
);
ntr
=
gridInqTrunc
(
gridID1
);
nlat
=
nfc2nlat
(
nfc
,
ntr
);
if
(
gridIDsp
!=
-
1
)
if
(
ntr
!=
gridInqTrunc
(
gridIDsp
)
)
gridIDsp
=
-
1
;
if
(
gridIDsp
==
-
1
)
{
nsp
=
(
ntr
+
1
)
*
(
ntr
+
2
);
gridIDsp
=
gridCreate
(
GRID_SPECTRAL
,
nsp
);
gridDefTrunc
(
gridIDsp
,
ntr
);
gridDefComplexPacking
(
gridIDsp
,
1
);
}
gridID2
=
gridIDsp
;
nlon
=
2
*
nlat
;
ntr
=
gridInqTrunc
(
gridID2
);
sptrans
=
sptrans_new
(
nlon
,
nlat
,
ntr
,
0
);
}
}
else
if
(
operatorID
==
SP2FC
)
{
if
(
gridIDsp
==
-
1
)
cdoWarning
(
"No spectral data found!"
);
gridID1
=
gridIDsp
;
if
(
gridID1
!=
-
1
)
{
ntr
=
gridInqTrunc
(
gridID1
);
nlat
=
ntr2nlat
(
ntr
);
if
(
gridIDfc
!=
-
1
)
{
if
(
ntr
!=
gridInqTrunc
(
gridIDfc
)
)
gridIDfc
=
-
1
;
}
if
(
gridIDfc
==
-
1
)
{
nfc
=
2
*
nlat
*
(
ntr
+
1
);
gridIDfc
=
gridCreate
(
GRID_FOURIER
,
nfc
);
gridDefTrunc
(
gridIDfc
,
ntr
);
}
gridID2
=
gridIDfc
;
nlon
=
2
*
nlat
;
sptrans
=
sptrans_new
(
nlon
,
nlat
,
ntr
,
0
);
}
}
else
if
(
operatorID
==
GP2FC
)
{
if
(
gridIDgp
==
-
1
)
cdoWarning
(
"No Gaussian grid data found!"
);
gridID1
=
gridIDgp
;
if
(
gridID1
!=
-
1
)
{
nlon
=
gridInqXsize
(
gridID1
);
nlat
=
gridInqYsize
(
gridID1
);
ntr
=
nlat2ntr
(
nlat
);
if
(
gridIDfc
!=
-
1
)
if
(
ntr
!=
gridInqTrunc
(
gridIDfc
)
)
gridIDfc
=
-
1
;
if
(
gridIDfc
==
-
1
)
{
nfc
=
2
*
nlat
*
(
ntr
+
1
);
gridIDfc
=
gridCreate
(
GRID_FOURIER
,
nfc
);
gridDefTrunc
(
gridIDfc
,
ntr
);
}
gridID2
=
gridIDfc
;
nfc
=
gridInqSize
(
gridID2
);
sptrans
=
sptrans_new
(
nlon
,
nlat
,
ntr
,
0
);
}
}
else
if
(
operatorID
==
FC2GP
)
{
if
(
gridIDfc
==
-
1
)
cdoWarning
(
"No fourier data found!"
);
gridID1
=
gridIDfc
;
if
(
gridID1
!=
-
1
)
{
nfc
=
gridInqSize
(
gridID1
);
ntr
=
gridInqTrunc
(
gridID1
);
nlat
=
nfc2nlat
(
nfc
,
ntr
);
if
(
gridIDgp
!=
-
1
)
{
if
(
nlat
!=
gridInqYsize
(
gridIDgp
)
)
gridIDgp
=
-
1
;
}
if
(
gridIDgp
==
-
1
)
{
char
gridname
[
20
];
sprintf
(
gridname
,
"t%dgrid"
,
ntr
);
gridIDgp
=
gridFromName
(
gridname
);
}
gridID2
=
gridIDgp
;
nlon
=
gridInqXsize
(
gridID2
);
nlat
=
gridInqYsize
(
gridID2
);
sptrans
=
sptrans_new
(
nlon
,
nlat
,
ntr
,
0
);
}
}
// printf("nfc %d, ntr %d, nlat %d, nlon %d\n", nfc, ntr, nlat, nlon);
nvars
=
vlistNvars
(
vlistID2
);
vars
=
(
int
*
)
malloc
(
nvars
*
sizeof
(
int
));
for
(
varID
=
0
;
varID
<
nvars
;
varID
++
)
{
if
(
gridID1
==
vlistInqVarGrid
(
vlistID1
,
varID
)
)
vars
[
varID
]
=
TRUE
;
else
vars
[
varID
]
=
FALSE
;
}
if
(
gridID1
!=
-
1
)
vlistChangeGrid
(
vlistID2
,
gridID1
,
gridID2
);
streamID2
=
streamOpenWrite
(
cdoStreamName
(
1
),
cdoFiletype
());
streamDefVlist
(
streamID2
,
vlistID2
);
gridsize
=
vlistGridsizeMax
(
vlistID1
);
array1
=
(
double
*
)
malloc
(
gridsize
*
sizeof
(
double
));
if
(
gridID2
!=
-
1
)
{
gridsize
=
gridInqSize
(
gridID2
);
array2
=
(
double
*
)
malloc
(
gridsize
*
sizeof
(
double
));
}
tsID
=
0
;
while
(
(
nrecs
=
streamInqTimestep
(
streamID1
,
tsID
))
)
{
taxisCopyTimestep
(
taxisID2
,
taxisID1
);
streamDefTimestep
(
streamID2
,
tsID
);
for
(
recID
=
0
;
recID
<
nrecs
;
recID
++
)
{
streamInqRecord
(
streamID1
,
&
varID
,
&
levelID
);
if
(
vars
[
varID
]
)
{
streamReadRecord
(
streamID1
,
array1
,
&
nmiss
);
if
(
nmiss
)
cdoAbort
(
"Missing values unsupported for spectral/fourier data!"
);
gridID1
=
vlistInqVarGrid
(
vlistID1
,
varID
);
if
(
operatorID
==
FC2SP
)
four2spec
(
sptrans
,
gridID1
,
array1
,
gridID2
,
array2
);
else
if
(
operatorID
==
SP2FC
)
spec2four
(
sptrans
,
gridID1
,
array1
,
gridID2
,
array2
);
else
if
(
operatorID
==
FC2GP
)
four2grid
(
sptrans
,
gridID1
,
array1
,
gridID2
,
array2
);
else
if
(
operatorID
==
GP2FC
)
grid2four
(
sptrans
,
gridID1
,
array1
,
gridID2
,
array2
);
streamDefRecord
(
streamID2
,
varID
,
levelID
);
streamWriteRecord
(
streamID2
,
array2
,
nmiss
);
}
else
{
streamDefRecord
(
streamID2
,
varID
,
levelID
);
if
(
lcopy
)
{
streamCopyRecord
(
streamID2
,
streamID1
);
}
else
{
streamReadRecord
(
streamID1
,
array1
,
&
nmiss
);
streamWriteRecord
(
streamID2
,
array1
,
nmiss
);
}
}
}
tsID
++
;
}
streamClose
(
streamID2
);
streamClose
(
streamID1
);
if
(
array2
)
free
(
array2
);
if
(
array1
)
free
(
array1
);
if
(
vars
)
free
(
vars
);
sptrans_delete
(
sptrans
);
cdoFinish
();
return
(
0
);
}
src/Makefile.am
View file @
19eb9b95
...
...
@@ -38,6 +38,7 @@ cdo_SOURCES += Arith.c \
Eofcoeff.c
\
Eofcoeff3d.c
\
Exprf.c
\
FC.c
\
Filedes.c
\
Fillmiss.c
\
Filter.c
\
...
...
src/Makefile.in
View file @
19eb9b95
...
...
@@ -71,7 +71,7 @@ am_cdo_OBJECTS = cdo-cdo.$(OBJEXT) cdo-Arith.$(OBJEXT) \
cdo-Enlarge.
$(OBJEXT)
cdo-Enlargegrid.
$(OBJEXT)
\
cdo-Ensstat.
$(OBJEXT)
cdo-Ensstat3.
$(OBJEXT)
\
cdo-Ensval.
$(OBJEXT)
cdo-Eofcoeff.
$(OBJEXT)
\
cdo-Eofcoeff3d.
$(OBJEXT)
cdo-Exprf.
$(OBJEXT)
\
cdo-Eofcoeff3d.
$(OBJEXT)
cdo-Exprf.
$(OBJEXT)
cdo-FC.
$(OBJEXT)
\
cdo-Filedes.
$(OBJEXT)
cdo-Fillmiss.
$(OBJEXT)
\
cdo-Filter.
$(OBJEXT)
cdo-Fldrms.
$(OBJEXT)
\
cdo-Fldstat.
$(OBJEXT)
cdo-Fldstat2.
$(OBJEXT)
\
...
...
@@ -345,18 +345,18 @@ cdo_SOURCES = cdo.c Arith.c Arithc.c Arithdays.c Arithlat.c CDItest.c \
Copy.c Deltime.c Derivepar.c Detrend.c Diff.c Duplicate.c
\
EOFs.c Eof3d.c EcaIndices.c Echam5ini.c Enlarge.c
\
Enlargegrid.c Ensstat.c Ensstat3.c Ensval.c Eofcoeff.c
\
Eofcoeff3d.c Exprf.c Filedes.c Fillmiss.c Filter.c
Fldrms.c
\
Fldstat.c Fldstat2.c Fourier.c Gather.c Gengrid.c
Gradsdes.c
\
Gridboxstat.c Gridcell.c Harmonic.c Hi.c
Histogram.c
\
IFS2ICON.c Importamsr.c Importbinary.c
Importcmsaf.c
\
Importobs.c Info.c Input.c Intgrid.c
Intgridtraj.c Intlevel.c
\
Intlevel3d.c Intntime.c Inttime.c
Intyear.c Invert.c
\
Invertlev.c Isosurface.c Log.c Maskbox.c
Mastrfu.c Math.c
\
Merge.c Mergegrid.c Mergetime.c Merstat.c
Monarith.c Mrotuv.c
\
Mrotuvb.c Ninfo.c Nmltest.c Output.c
Outputgmt.c Pinfo.c
\
Pressure.c Regres.c Remap.c Remapeta.c
Replace.c
\
Replacevalues.c Rotuv.c Runpctl.c Runstat.c
Scatter.c
\
Seascount.c Seaspctl.c Seasstat.c Selbox.c Select.c
\
Eofcoeff3d.c Exprf.c
FC.c
Filedes.c Fillmiss.c Filter.c
\
Fldrms.c
Fldstat.c Fldstat2.c Fourier.c Gather.c Gengrid.c
\
Gradsdes.c
Gridboxstat.c Gridcell.c Harmonic.c Hi.c
\
Histogram.c
IFS2ICON.c Importamsr.c Importbinary.c
\
Importcmsaf.c
Importobs.c Info.c Input.c Intgrid.c
\
Intgridtraj.c Intlevel.c
Intlevel3d.c Intntime.c Inttime.c
\
Intyear.c Invert.c
Invertlev.c Isosurface.c Log.c Maskbox.c
\
Mastrfu.c Math.c
Merge.c Mergegrid.c Mergetime.c Merstat.c
\
Monarith.c Mrotuv.c
Mrotuvb.c Ninfo.c Nmltest.c Output.c
\
Outputgmt.c Pinfo.c
Pressure.c Regres.c Remap.c Remapeta.c
\
Replace.c
Replacevalues.c Rotuv.c Runpctl.c Runstat.c
\
Scatter.c
Seascount.c Seaspctl.c Seasstat.c Selbox.c Select.c
\
Seloperator.c Selrec.c Seltime.c Selvar.c Set.c Setbox.c
\
Setgatt.c Setgrid.c Sethalo.c Setmiss.c Setrcaname.c Settime.c
\
Setzaxis.c Showinfo.c Sinfo.c Smooth9.c Sort.c Sorttimestamp.c
\
...
...
@@ -549,6 +549,7 @@ distclean-compile:
@AMDEP_TRUE@@am__include@
@am__quote@./$(DEPDIR)/cdo-Eofcoeff.Po@am__quote@
@AMDEP_TRUE@@am__include@
@am__quote@./$(DEPDIR)/cdo-Eofcoeff3d.Po@am__quote@
@AMDEP_TRUE@@am__include@
@am__quote@./$(DEPDIR)/cdo-Exprf.Po@am__quote@
@AMDEP_TRUE@@am__include@
@am__quote@./$(DEPDIR)/cdo-FC.Po@am__quote@
@AMDEP_TRUE@@am__include@
@am__quote@./$(DEPDIR)/cdo-Filedes.Po@am__quote@
@AMDEP_TRUE@@am__include@
@am__quote@./$(DEPDIR)/cdo-Fillmiss.Po@am__quote@
@AMDEP_TRUE@@am__include@
@am__quote@./$(DEPDIR)/cdo-Filter.Po@am__quote@
...
...
@@ -1261,6 +1262,20 @@ cdo-Exprf.obj: Exprf.c
@AMDEP_TRUE@@am__fastdepCC_FALSE@
DEPDIR
=
$(DEPDIR)
$(CCDEPMODE)
$(depcomp)
@AMDEPBACKSLASH@
@am__fastdepCC_FALSE@
$(CC)
$(DEFS)
$(DEFAULT_INCLUDES)
$(INCLUDES)
$(cdo_CPPFLAGS)
$(CPPFLAGS)
$(AM_CFLAGS)
$(CFLAGS)
-c
-o
cdo-Exprf.obj
`if
test
-f
'Exprf.c'
;
then
$(CYGPATH_W)
'Exprf.c'
;
else
$(CYGPATH_W)
'$(srcdir)/Exprf.c'
;
fi`
cdo-FC.o
:
FC.c
@am__fastdepCC_TRUE@
$(CC)
$(DEFS)
$(DEFAULT_INCLUDES)
$(INCLUDES)
$(cdo_CPPFLAGS)
$(CPPFLAGS)
$(AM_CFLAGS)
$(CFLAGS)
-MT
cdo-FC.o
-MD
-MP
-MF
$(DEPDIR)/cdo-FC.Tpo
-c
-o
cdo-FC.o
`test
-f
'FC.c'
||
echo
'$(srcdir)/'
`FC.c
@am__fastdepCC_TRUE@
$(am__mv)
$(DEPDIR)/cdo-FC.Tpo
$(DEPDIR)/cdo-FC.Po
@AMDEP_TRUE@@am__fastdepCC_FALSE@
source
=
'FC.c'
object
=
'cdo-FC.o'
libtool
=
no @AMDEPBACKSLASH@
@AMDEP_TRUE@@am__fastdepCC_FALSE@
DEPDIR
=
$(DEPDIR)
$(CCDEPMODE)
$(depcomp)
@AMDEPBACKSLASH@
@am__fastdepCC_FALSE@
$(CC)
$(DEFS)
$(DEFAULT_INCLUDES)
$(INCLUDES)
$(cdo_CPPFLAGS)
$(CPPFLAGS)
$(AM_CFLAGS)
$(CFLAGS)
-c
-o
cdo-FC.o
`test
-f
'FC.c'
||
echo
'$(srcdir)/'
`FC.c
cdo-FC.obj
:
FC.c
@am__fastdepCC_TRUE@
$(CC)
$(DEFS)
$(DEFAULT_INCLUDES)
$(INCLUDES)
$(cdo_CPPFLAGS)
$(CPPFLAGS)
$(AM_CFLAGS)
$(CFLAGS)
-MT
cdo-FC.obj
-MD
-MP
-MF
$(DEPDIR)/cdo-FC.Tpo
-c
-o
cdo-FC.obj
`if
test
-f
'FC.c'
;
then
$(CYGPATH_W)
'FC.c'
;
else
$(CYGPATH_W)
'$(srcdir)/FC.c'
;
fi`
@am__fastdepCC_TRUE@
$(am__mv)
$(DEPDIR)/cdo-FC.Tpo
$(DEPDIR)/cdo-FC.Po
@AMDEP_TRUE@@am__fastdepCC_FALSE@
source
=
'FC.c'
object
=
'cdo-FC.obj'
libtool
=
no @AMDEPBACKSLASH@
@AMDEP_TRUE@@am__fastdepCC_FALSE@
DEPDIR
=
$(DEPDIR)
$(CCDEPMODE)
$(depcomp)
@AMDEPBACKSLASH@
@am__fastdepCC_FALSE@
$(CC)
$(DEFS)
$(DEFAULT_INCLUDES)
$(INCLUDES)
$(cdo_CPPFLAGS)
$(CPPFLAGS)
$(AM_CFLAGS)
$(CFLAGS)
-c
-o
cdo-FC.obj
`if
test
-f
'FC.c'
;
then
$(CYGPATH_W)
'FC.c'
;
else
$(CYGPATH_W)
'$(srcdir)/FC.c'
;
fi`
cdo-Filedes.o
:
Filedes.c
@am__fastdepCC_TRUE@
$(CC)
$(DEFS)
$(DEFAULT_INCLUDES)
$(INCLUDES)
$(cdo_CPPFLAGS)
$(CPPFLAGS)
$(AM_CFLAGS)
$(CFLAGS)
-MT
cdo-Filedes.o
-MD
-MP
-MF
$(DEPDIR)/cdo-Filedes.Tpo
-c
-o
cdo-Filedes.o
`test
-f
'Filedes.c'
||
echo
'$(srcdir)/'
`Filedes.c
@am__fastdepCC_TRUE@
$(am__mv)
$(DEPDIR)/cdo-Filedes.Tpo
$(DEPDIR)/cdo-Filedes.Po
...
...
src/Spectral.c
View file @
19eb9b95
...
...
@@ -260,7 +260,7 @@ void *Spectral(void *argument)
if
(
vars
[
varID
]
)
{
streamReadRecord
(
streamID1
,
array1
,
&
nmiss
);
if
(
nmiss
)
cdoAbort
(
"
m
issing values unsupported for spectral data!"
);
if
(
nmiss
)
cdoAbort
(
"
M
issing values unsupported for spectral data!"
);
gridID1
=
vlistInqVarGrid
(
vlistID1
,
varID
);
if
(
operatorID
==
GP2SP
||
operatorID
==
GP2SPL
)
...
...
src/cdo_int.h
View file @
19eb9b95
...
...
@@ -124,6 +124,8 @@ int zaxis2ltype(int zaxisID);
int
ztype2ltype
(
int
zaxistype
);
int
ltype2ztype
(
int
ltype
);
int
nfc2nlat
(
int
nfc
,
int
ntr
);
int
nlat2ntr
(
int
nlat
);
int
nlat2ntr_linear
(
int
nlat
);
int
ntr2nlat
(
int
ntr
);
...
...
src/griddes.c
View file @
19eb9b95
...
...
@@ -1357,6 +1357,17 @@ int gridFromPingo(FILE *gfp, const char *dname)
}
int
nfc2nlat
(
int
nfc
,
int
ntr
)
{
int
nlat
;
nlat
=
nfc
/
(
ntr
+
1
);
nlat
/=
2
;
return
(
nlat
);
}
int
nlat2ntr
(
int
nlat
)
{
int
ntr
;
...
...
src/modules.c
View file @
19eb9b95
...
...
@@ -72,6 +72,7 @@ void *Eofcoeff3d(void *argument);
void
*
EOFs
(
void
*
argument
);
void
*
EOF3d
(
void
*
argument
);
void
*
Expr
(
void
*
argument
);
void
*
FC
(
void
*
argument
);
void
*
Filedes
(
void
*
argument
);
void
*
Fillmiss
(
void
*
argument
);
void
*
Filter
(
void
*
argument
);
...
...
@@ -283,6 +284,7 @@ void *Wct(void *argument);
#define EOFsOperators {"eof", "eofspatial", "eoftime"}
#define EOF3dOperators {"eof3d","eof3dspatial","eof3dtime"}
#define ExprOperators {"expr", "exprf", "aexpr", "aexprf"}
#define FCOperators {"fc2sp", "sp2fc", "fc2gp", "gp2fc"}
#define FiledesOperators {"filedes", "griddes", "griddes2", "zaxisdes", "vct", "vct2", "pardes", \
"vlist", "partab", "partab2"}
#define FillmissOperators {"fillmiss"}
...
...
@@ -529,6 +531,7 @@ static modules_t Modules[] =
{
EOFs
,
EOFsHelp
,
EOFsOperators
,
CDI_REAL
,
1
,
2
},
{
EOF3d
,
NULL
,
EOF3dOperators
,
CDI_REAL
,
1
,
2
},
{
Expr
,
ExprHelp
,
ExprOperators
,
CDI_REAL
,
1
,
1
},
{
FC
,
NULL
,
FCOperators
,
CDI_REAL
,
1
,
1
},
{
Filedes
,
FiledesHelp
,
FiledesOperators
,
CDI_BOTH
,
1
,
0
},
{
Fillmiss
,
NULL
,
FillmissOperators
,
CDI_REAL
,
1
,
1
},
{
Filter
,
FilterHelp
,
FilterOperators
,
CDI_REAL
,
1
,
1
},
...
...
src/printinfo.h
View file @
19eb9b95
...
...
@@ -249,14 +249,13 @@ void printGridInfo(int vlistID)
}
else
if
(
gridtype
==
GRID_SPECTRAL
)
{
fprintf
(
stdout
,
"size : dim = %d truncation = %d nsp = %d
\n
"
,
gridsize
,
trunc
,
gridsize
/
2
);
fprintf
(
stdout
,
"size : dim = %d nsp = %d truncation = %d
\n
"
,
gridsize
,
gridsize
/
2
,
trunc
);
fprintf
(
stdout
,
"%*s"
,
nbyte0
,
""
);
fprintf
(
stdout
,
" complexPacking = %d
\n
"
,
gridInqComplexPacking
(
gridID
));
}
else
if
(
gridtype
==
GRID_FOURIER
)
{
fprintf
(
stdout
,
"size : dim = %d nfc = %d
\n
"
,
gridsize
,
gridsize
/
2
);
fprintf
(
stdout
,
"size : dim = %d nfc =
%d truncation =
%d
\n
"
,
gridsize
,
gridsize
/
2
,
trunc
);
}
else
if
(
gridtype
==
GRID_GME
)
{
...
...
src/specspace.c
View file @
19eb9b95
...
...
@@ -2,8 +2,9 @@
#include
<math.h>
#include
<string.h>
#include
"cdo.h"
#include
<cdi.h>
#include
"cdo.h"
#include
"cdo_int.h"
#ifndef _DMEMORY_H
# include "dmemory.h"
#endif
...
...
@@ -220,6 +221,73 @@ void spec2grid(SPTRANS *sptrans, int gridIDin, double *arrayIn, int gridIDout, d
}
void
four2spec
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
)
{
int
ntr
,
nlat
,
nfc
;
int
nlev
=
1
;
int
waves
;
ntr
=
gridInqTrunc
(
gridIDout
);
nlat
=
sptrans
->
nlat
;
waves
=
ntr
+
1
;
nfc
=
waves
*
2
;
fc2sp
(
arrayIn
,
arrayOut
,
sptrans
->
pold
,
nlev
,
nlat
,
nfc
,
ntr
);
}
void
spec2four
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
)
{
int
ntr
,
nlat
,
nfc
;
int
nlev
=
1
;
int
waves
;
ntr
=
gridInqTrunc
(
gridIDin
);
nfc
=
gridInqSize
(
gridIDout
);
nlat
=
nfc2nlat
(
nfc
,
ntr
);
waves
=
ntr
+
1
;
nfc
=
waves
*
2
;
sp2fc
(
arrayIn
,
arrayOut
,
sptrans
->
poli
,
nlev
,
nlat
,
nfc
,
ntr
);
}
void
four2grid
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
)
{
int
ntr
,
nlat
,
nlon
,
nfc
;
int
nlev
=
1
;
int
waves
;
ntr
=
gridInqTrunc
(
gridIDin
);
nlon
=
gridInqXsize
(
gridIDout
);
nlat
=
gridInqYsize
(
gridIDout
);
waves
=
ntr
+
1
;
nfc
=
waves
*
2
;
fc2gp
(
sptrans
->
trig
,
sptrans
->
ifax
,
arrayIn
,
arrayOut
,
nlat
,
nlon
,
nlev
,
nfc
);
}
void
grid2four
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
)
{
int
nlat
,
nlon
,
nfc
,
ntr
;
int
nlev
=
1
;
int
waves
;
ntr
=
gridInqTrunc
(
gridIDout
);
nlon
=
gridInqXsize
(
gridIDin
);
nlat
=
gridInqYsize
(
gridIDin
);
waves
=
ntr
+
1
;
nfc
=
waves
*
2
;
gp2fc
(
sptrans
->
trig
,
sptrans
->
ifax
,
arrayIn
,
arrayOut
,
nlat
,
nlon
,
nlev
,
nfc
);
}
void
spec2spec
(
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
)
{
int
ntrIn
,
ntrOut
;
...
...
@@ -244,7 +312,7 @@ void speccut(int gridIDin, double *arrayIn, double *arrayOut, int *waves)
SPTRANS
*
sptrans_new
(
int
nlon
,
int
nlat
,
int
ntr
,
int
flag
)
{
SPTRANS
*
sptrans
;
int
dim
sp
;
int
n
sp
;
sptrans
=
(
SPTRANS
*
)
malloc
(
sizeof
(
SPTRANS
));
...
...
@@ -252,8 +320,8 @@ SPTRANS *sptrans_new(int nlon, int nlat, int ntr, int flag)
sptrans
->
nlat
=
nlat
;
sptrans
->
ntr
=
ntr
;
dim
sp
=
(
ntr
+
1
)
*
(
ntr
+
2
);
sptrans
->
poldim
=
dim
sp
/
2
*
nlat
;
n
sp
=
(
ntr
+
1
)
*
(
ntr
+
2
);
sptrans
->
poldim
=
n
sp
/
2
*
nlat
;
sptrans
->
trig
=
(
double
*
)
malloc
(
nlon
*
sizeof
(
double
));
fft_set
(
sptrans
->
trig
,
sptrans
->
ifax
,
nlon
);
...
...
src/specspace.h
View file @
19eb9b95
...
...
@@ -41,6 +41,10 @@ void trans_dv2uv(SPTRANS *sptrans, DVTRANS *dvtrans, int nlev,
void
grid2spec
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
);
void
spec2grid
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
);
void
four2spec
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
);
void
spec2four
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
);
void
four2grid
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
);
void
grid2four
(
SPTRANS
*
sptrans
,
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
);
void
spec2spec
(
int
gridIDin
,
double
*
arrayIn
,
int
gridIDout
,
double
*
arrayOut
);
void
speccut
(
int
gridIDin
,
double
*
arrayIn
,
double
*
arrayOut
,
int
*
waves
);
...
...
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment