#include #include #include #include #include #include #include #include #include #include #include #include "cctk.h" #include "cctk_Parameters.h" #include "AMRwriter.hh" #include "AmrGridReader.hh" #ifdef HDF4 # include "HDFIO.hh" #endif #ifdef HDF5 # include "H5IO.hh" #endif #include "IEEEIO.hh" #include "IO.hh" // Hack to stop FlexIO data type clash with LAM MPI #undef BYTE #undef CHAR #include "CactusBase/IOUtil/src/ioGH.h" #include "bbox.hh" #include "data.hh" #include "gdata.hh" #include "ggf.hh" #include "vect.hh" #include "carpet.hh" #include "ioflexio.hh" extern "C" { static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/CarpetAttic/CarpetIOFlexIOCheckpoint/src/ioflexio.cc,v 1.13 2003/12/10 14:49:10 cott Exp $"; CCTK_FILEVERSION(Carpet_CarpetIOFlexIO_ioflexio_cc); } namespace CarpetIOFlexIO { using namespace std; using namespace Carpet; using namespace CarpetIOFlexIOUtil; // Variable definitions // int GHExtension; int IOMethod; vector do_truncate; vector > last_output; static const char* GetStringParameter (const char* const parametername, const char* const fallback); static int GetIntParameter (const char* const parametername, int fallback); static bool CheckForVariable (const cGH* const cgh, const char* const varlist, const int vindex); static void SetFlag (int index, const char* optstring, void* arg); int CarpetIOFlexIO_Startup () { CCTK_RegisterBanner ("AMR 3D FlexIO I/O provided by CarpetIOFlexIO"); // GHExtension = CCTK_RegisterGHExtension("CarpetIOFlexIO"); CCTK_RegisterGHExtensionSetupGH (CCTK_RegisterGHExtension("CarpetIOFlexIO"),SetupGH); IOMethod = CCTK_RegisterIOMethod ("CarpetIOFlexIO"); CCTK_RegisterIOMethodOutputGH (IOMethod, OutputGH); CCTK_RegisterIOMethodOutputVarAs (IOMethod, OutputVarAs); CCTK_RegisterIOMethodTimeToOutput (IOMethod, TimeToOutput); CCTK_RegisterIOMethodTriggerOutput (IOMethod, TriggerOutput); return 0; } void* SetupGH (tFleshConfig* const fc, const int convLevel, cGH* const cgh) { DECLARE_CCTK_PARAMETERS; CarpetIOFlexIOGH* myGH; CCTK_INT i; // Truncate all files if this is not a restart do_truncate.resize(CCTK_NumVars(), true); // No iterations have yet been output last_output.resize(maxreflevels); for (int rl=0; rlout_last = (int *) malloc (numvars * sizeof (int)); myGH->requests = (ioRequest **) calloc (numvars, sizeof (ioRequest *)); myGH->cp_filename_list = (char **) calloc (abs (checkpoint_keep), sizeof (char *)); myGH->cp_filename_index = 0; myGH->out_vars = strdup (""); myGH->out_every_default = out_every - 1; for (i = 0; i < numvars; i++) { myGH->out_last [i] = -1; } myGH->open_output_files = NULL; return (myGH); return 0; } int OutputGH (const cGH* const cgh) { for (int vindex=0; vindexvindex; const int group = CCTK_GroupIndexFromVarI (varindex); const int n0 = CCTK_FirstVarIndexI(group); assert (n0>=0 && n0=0 && var0)); // if(grouptype == CCTK_SCALAR || grouptype == CCTK_ARRAY) return 0; if (CCTK_MyProc(cgh)==0) { // Set datatype //fprintf(stderr,"%d\n",CCTK_VarTypeI(varindex)); /* should be obsolete now assert (CCTK_VarTypeI(varindex) == CCTK_VARIABLE_REAL8 || (sizeof(CCTK_REAL) == sizeof(CCTK_REAL8) && CCTK_VarTypeI(varindex) == CCTK_VARIABLE_REAL)); */ amrwriter->setType (FlexIODataType(CCTK_VarTypeI(varindex))); // Set name information // char *name = CCTK_FullName (varindex); //IOwriteAttribute(amrwriter->file,"name",IObase::Char,strlen(name)+1,name); //free(name); int gpdim = CCTK_GroupDimI(group); // need gpdim=1 if scalar (flexio wants this) if(gpdim == 0) if(grouptype == CCTK_SCALAR) gpdim = 1; else CCTK_WARN(0,"Non-scalar variable with dimension 0!!!"); // Set coordinate information CCTK_REAL lower[dim], upper[dim]; double origin[dim], delta[dim], timestep; for (int d=0; dcctk_delta_space[d]; } timestep = cgh->cctk_delta_time; amrwriter->setTopLevelParameters (gpdim, origin, delta, timestep, maxreflevels); // Set refinement information int interlevel_timerefinement; int interlevel_spacerefinement[dim]; int initial_gridplacementrefinement[dim]; interlevel_timerefinement = hh->reffact; for (int d=0; dreffact; initial_gridplacementrefinement[d] = 1; } amrwriter->setRefinement (interlevel_timerefinement, interlevel_spacerefinement, initial_gridplacementrefinement); // Set level amrwriter->setLevel (reflevel); // Set current time amrwriter->setTime (cgh->cctk_iteration); } // Traverse all components on this refinement and multigrid level BEGIN_COMPONENT_LOOP(cgh, grouptype) { const ggf* ff = 0; #if 0 if (grouptype == CCTK_ARRAY){ // this is a DIRTY hack to fix problems caused by the fact that I am to lazy to write a more // general output routine... if(reflevel !=0) return 0; } else #endif //if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "GF reflevel: %d component: %d grouptype: %d",reflevel,component,grouptype); assert (var < (int)arrdata[group].data.size()); ff = (ggf*)arrdata[group].data[var]; // CCTK_VInfo (CCTK_THORNSTRING, "bogus reflevel,component,mglevel %d,%d,%d",reflevel,component,mglevel); const gdata* const data = (*ff) (tl, reflevel, component, mglevel); // get some more group information cGroupDynamicData gdyndata; int ierr = CCTK_GroupDynamicData(cgh,group,&gdyndata); assert(ierr==0); // Make temporary copy on processor 0 bbox ext = data->extent(); vect lo = ext.lower(); vect hi = ext.upper(); vect str = ext.stride(); // Ignore ghost zones if desired #warning "need to check size of gdyndata.bbox" for (int d=0; d(lo,hi,str); gdata* const tmp = data->make_typed (varindex); tmp->allocate (ext, 0); for (comm_state state; !state.done(); state.step()) { tmp->copy_from (state, data, ext); } // Write data if (CCTK_MyProc(cgh)==0) { int origin[dim], dims[dim]; for (int d=0; dwrite (origin, dims, (void*)tmp->storage()); // dump attributes DumpCommonAttributes(cgh,writer,request); // char *name = CCTK_FullName (varindex); // writer->writeAttribute("name",IObase::Char,strlen(name)+1,name); //free(name); } // Delete temporary copy delete tmp; } END_COMPONENT_LOOP; return 0; } int WriteGS (const cGH* const cgh, IObase* writer, ioRequest* request) { #warning This function should be obsolete by now!!! // writes out a grid scalar DECLARE_CCTK_PARAMETERS; const int timelevel = request->timelevel; const int varindex = request->vindex; const int group = CCTK_GroupIndexFromVarI (varindex); const int n0 = CCTK_FirstVarIndexI(group); assert (n0>=0 && n0=0 && var0)); int myproc = CCTK_MyProc (cgh); int nprocs = CCTK_nProcs (cgh); char* fullname = CCTK_FullName (varindex); int datatype = CCTK_VarTypeI(varindex); int datatypesize = CCTK_VarTypeSize(datatype); char* buffer = (char*) calloc (nprocs, datatypesize); memcpy (buffer + myproc*datatypesize, CCTK_VarDataPtrI (cgh,timelevel,varindex), datatypesize); if (nprocs > 1) { int i = CCTK_ReductionHandle ("sum"); if (i >= 0) { i = CCTK_ReduceArray (cgh, -1, i, nprocs, datatype, buffer, 1, 1, datatype, nprocs, buffer); } if (i < 0) { CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING, "WriteGS: Cannot check whether values on differentprocessors are the same for grid scalar '%s'", fullname); // copy this processor's value to the start of buffer memcpy (buffer, buffer + myproc*datatypesize, datatypesize); } else { int retval = 0; for (i = 1; i < nprocs; i++) { retval |= memcmp (buffer, buffer + i*datatypesize, datatypesize); } if (retval) { CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING, "WriteGS: value of grid scalar variable '%s' (timelevel %d)" " differs between processors, only value from processor 0 " "will be written", fullname, timelevel); } } } if (myproc==0) { int dim = 1; // Traverse all components on this refinement and multigrid level BEGIN_COMPONENT_LOOP(cgh, grouptype) { // actually, looping makes no sense here, since a scalar must be the // same on all components. in fact, the loop is not being // executed for scalars; see macro definition. if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "SCALAR reflevel,component,mglevel %d,%d,%d",reflevel,component,mglevel); writer->write(FlexIODataType(CCTK_VarTypeI(varindex)),1,&dim,buffer); /* scalars have size 0 */ request->hsize[0] = 0; DumpCommonAttributes (cgh,writer,request); } END_COMPONENT_LOOP; } return 0; } int OutputVarAs (const cGH* const cgh, const char* const varname, const char* const alias) { DECLARE_CCTK_PARAMETERS; const int n = CCTK_VarIndex(varname); assert (n>=0 && n=0 && group<(int)Carpet::arrdata.size()); const int n0 = CCTK_FirstVarIndexI(group); assert (n0>=0 && n0=0 && var0) return 0; int first_vindex = CCTK_FirstVarIndexI (group); /* get the default I/O request for this group */ ioRequest* request = IOUtil_DefaultIORequest (cgh, first_vindex, 1); // Get grid hierarchy extentsion from IOUtil const ioGH * const iogh = (const ioGH *)CCTK_GHExtension (cgh, "IO"); assert (iogh); // Create the output directory const char* myoutdir = GetStringParameter("out3D_dir", out_dir); if (CCTK_MyProc(cgh)==0) { CCTK_CreateDirectory (0755, myoutdir); } // Invent a file name const char* extension = 0; if (CCTK_Equals(out3D_format, "IEEE")) { extension = ".raw"; #ifdef HDF4 } else if (CCTK_Equals(out3D_format, "HDF4")) { extension = ".hdf"; #endif #ifdef HDF5 } else if (CCTK_Equals(out3D_format, "HDF5")) { extension = ".h5"; #endif } else { assert (0); } extension = GetStringParameter ("out3D_extension", extension); ostringstream filenamebuf; filenamebuf << myoutdir << "/" << alias << extension; string filenamestr = filenamebuf.str(); const char * const filename = filenamestr.c_str(); IObase* writer = 0; AMRwriter* amrwriter = 0; // Write the file only on the root processor if (CCTK_MyProc(cgh)==0) { // If this is the first time, then create and truncate the file if (do_truncate[n]) { struct stat fileinfo; if (! iogh->recovered || stat(filename, &fileinfo)!=0) { writer = 0; if (CCTK_Equals(out3D_format, "IEEE")) { writer = new IEEEIO(filename, IObase::Create); #ifdef HDF4 } else if (CCTK_Equals(out3D_format, "HDF4")) { writer = new HDFIO(filename, IObase::Create); #endif #ifdef HDF5 } else if (CCTK_Equals(out3D_format, "HDF5")) { writer = new H5IO(filename, IObase::Create); #endif } else { assert (0); } delete writer; writer = 0; } } // Open the file if (CCTK_Equals(out3D_format, "IEEE")) { writer = new IEEEIO(filename, IObase::Append); #ifdef HDF4 } else if (CCTK_Equals(out3D_format, "HDF4")) { writer = new HDFIO(filename, IObase::Append); #endif #ifdef HDF5 } else if (CCTK_Equals(out3D_format, "HDF5")) { writer = new H5IO(filename, IObase::Append); #endif } else { assert (0); } assert (writer->isValid()); amrwriter = new AMRwriter(*writer); } WriteGF(cgh,writer,amrwriter,request); // Close the file if (CCTK_MyProc(cgh)==0) { delete amrwriter; amrwriter = 0; delete writer; writer = 0; } // Don't truncate again do_truncate[n] = false; return 0; } int TimeToOutput (const cGH* const cgh, const int vindex) { DECLARE_CCTK_PARAMETERS; assert (vindex>=0 && vindex<(int)last_output[reflevel].size()); const int myoutevery = GetIntParameter("out3D_every", out_every); if (myoutevery < 0) { // Nothing should be output at all return 0; } if (cgh->cctk_iteration % myoutevery != 0) { // Nothing should be output during this iteration return 0; } if (! CheckForVariable(cgh, GetStringParameter("out3D_vars",""), vindex)) { // This variable should not be output return 0; } if (last_output[reflevel][vindex] == cgh->cctk_iteration) { // Has already been output during this iteration char* varname = CCTK_FullName(vindex); CCTK_VWarn (5, __LINE__, __FILE__, CCTK_THORNSTRING, "Skipping output for variable \"%s\", because this variable " "has already been output during the current iteration -- " "probably via a trigger during the analysis stage", varname); free (varname); return 0; } assert (last_output[reflevel][vindex] < cgh->cctk_iteration); // Should be output during this iteration return 1; } int TriggerOutput (const cGH* const cgh, const int vindex) { assert (vindex>=0 && vindexcctk_iteration; return retval; } int InputGH (const cGH* const cgh) { int retval = 0; for (int vindex=0; vindex=0 && n=0 && group<(int)Carpet::arrdata.size()); const int n0 = CCTK_FirstVarIndexI(group); assert (n0>=0 && n0=0 && varisValid()) { CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING, "Could not open file \"%s\" for reading", filename); } assert (reader->isValid()); if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Reading AMR info"); amrreader = new AmrGridReader(*reader); // Read information about dataset if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Reading dataset info"); IObase::DataType numbertype; reader->readInfo (numbertype, rank, dims); nbytes = IObase::nBytes(numbertype,rank,dims); if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "type=%d rank=%d dims=[%d,%d,%d] nbytes=%d", (int)numbertype, rank, dims[0], dims[1], dims[2], nbytes); // Check rank assert (rank==gpdim); // Check datatype // TODO: Check datatype correctly assert (CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL8 || (sizeof(CCTK_REAL) == sizeof(CCTK_REAL8) && CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL)); // TODO: check grid spacing // Number of datasets if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Reading number of datasets"); ndatasets = reader->nDatasets(); if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "ndatasets=%d", ndatasets); assert (ndatasets>=0); } // Broadcast rank, dimensions, and nbytes MPI_Bcast (&rank, 1, MPI_INT, 0, dist::comm); assert (rank>=1); MPI_Bcast (&dims, rank, MPI_INT, 0, dist::comm); for (int d=0; d=0); MPI_Bcast (&nbytes, 1, MPI_INT, 0, dist::comm); assert (nbytes>=0); // Broadcast number of datasets MPI_Bcast (&ndatasets, 1, MPI_INT, 0, dist::comm); assert (ndatasets>=0); // Read all datasets // TODO: read only some datasets for (int dataset=0; datasetgetGrid(dataset); assert (amrgrid!=0); assert (amrgrid->data!=0); // If iorigin attribute is absent, assume file has unigrid // data. Initialize iorigin to 0. IObase::DataType atype; int alength; if (reader->readAttributeInfo("iorigin", atype, alength) < 0) { for (int d=0; diorigin[d] = 0; } } for (int d=0; diorigin[d]; amr_dims[d] = amrgrid->dims[d]; } for (int d=gpdim; d* ff = 0; assert (var < (int)arrdata[group].data.size()); ff = (ggf*)arrdata[group].data[var]; gdata* const data = (*ff) (tl, rl, component, mglevel); // Create temporary data storage on processor 0 const vect str = vect(reflevelfact); const vect lb = vect(amr_origin) * str; const vect ub = lb + (vect(amr_dims) - 1) * str; const bbox ext(lb,ub,str); gdata* const tmp = data->make_typed (n); if (CCTK_MyProc(cgh)==0) { tmp->allocate (ext, 0, amrgrid->data); } else { tmp->allocate (ext, 0, 0); } // Copy into grid function for (comm_state state; !state.done(); state.step()) { data->copy_from (state, tmp, ext); } // Delete temporary copy delete tmp; } END_COMPONENT_LOOP; if (CCTK_MyProc(cgh)==0) { free (amrgrid->data); free (amrgrid); amrgrid = 0; } } // loop over datasets // Close the file if (CCTK_MyProc(cgh)==0) { if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Deleting AMR info"); delete amrreader; amrreader = 0; if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Closing file"); delete reader; reader = 0; } return 0; } int CarpetIOFlexIO_ReadData (CCTK_ARGUMENTS) { DECLARE_CCTK_ARGUMENTS; return InputGH(cctkGH); } const char* GetStringParameter (const char* const parametername, const char* const fallback) { if (CCTK_ParameterQueryTimesSet (parametername, CCTK_THORNSTRING) > 0) { int ptype; const char* const* const ppval = (const char* const*)CCTK_ParameterGet (parametername, CCTK_THORNSTRING, &ptype); assert (ppval); const char* const pval = *ppval; assert (ptype == PARAMETER_STRING); return pval; } return fallback; } int GetIntParameter (const char* const parametername, int fallback) { if (CCTK_ParameterQueryTimesSet (parametername, CCTK_THORNSTRING) > 0) { int ptype; const int* const ppval = (const int*)CCTK_ParameterGet (parametername, CCTK_THORNSTRING, &ptype); assert (ppval); const int pval = *ppval; assert (ptype == PARAMETER_INT); return pval; } return fallback; } bool CheckForVariable (const cGH* const cgh, const char* const varlist, const int vindex) { const int numvars = CCTK_NumVars(); assert (vindex>=0 && vindex flags(numvars); CCTK_TraverseString (varlist, SetFlag, &flags, CCTK_GROUP_OR_VAR); return flags[vindex]; } void SetFlag (int index, const char* optstring, void* arg) { vector& flags = *(vector*)arg; flags[index] = true; } } // namespace CarpetIOFlexIO