source: XIOS3/trunk/src/io/nc4_data_output.cpp @ 2386

Last change on this file since 2386 was 2386, checked in by jderouillat, 2 years ago

Set the code structure to compute the hash value of an element based on its attributs, use for now before writing an element in a file

  • Property copyright set to
    Software name : XIOS (Xml I/O Server)
    http://forge.ipsl.jussieu.fr/ioserver
    Creation date : January 2009
    Licence : CeCCIL version2
    see license file in root directory : Licence_CeCILL_V2-en.txt
    or http://www.cecill.info/licences/Licence_CeCILL_V2-en.html
    Holder : CEA/LSCE (Laboratoire des Sciences du CLimat et de l'Environnement)
    CNRS/IPSL (Institut Pierre Simon Laplace)
    Project Manager : Yann Meurdesoif
    yann.meurdesoif@cea.fr
  • Property svn:executable set to *
File size: 124.9 KB
Line 
1#include "nc4_data_output.hpp"
2
3#include "attribute_template.hpp"
4#include "group_template.hpp"
5
6#include "file.hpp"
7#include "calendar.hpp"
8#include "context.hpp"
9#include "context_server.hpp"
10#include "netCdfException.hpp"
11#include "exception.hpp"
12#include "timer.hpp"
13#include "uuid.hpp"
14namespace xios
15{
16      /// ////////////////////// Dfinitions ////////////////////// ///
17      CNc4DataOutput::CNc4DataOutput
18         (CFile* file, const StdString & filename, bool exist)
19            : SuperClass()
20            , SuperClassWriter(filename, exist)
21            , filename(filename)
22            , file(file),hasTimeInstant(false),hasTimeCentered(false), timeCounterType(none), relElements_()
23      {
24        SuperClass::type = MULTI_FILE;
25        compressionLevel= file->compression_level.isEmpty() ? 0 :file->compression_level ;
26      }
27
28      CNc4DataOutput::CNc4DataOutput
29         (CFile* file, const StdString & filename, bool exist, bool useClassicFormat, bool useCFConvention,
30          MPI_Comm comm_file, bool multifile, bool isCollective, const StdString& timeCounterName)
31            : SuperClass()
32            , SuperClassWriter(filename, exist, useClassicFormat, useCFConvention, &comm_file, multifile, timeCounterName)
33            , comm_file(comm_file)
34            , filename(filename)
35            , isCollective(isCollective)
36            , file(file),hasTimeInstant(false),hasTimeCentered(false), timeCounterType(none), relElements_()
37      {
38        SuperClass::type = (multifile) ? MULTI_FILE : ONE_FILE;
39        if (file==NULL) compressionLevel = 0 ;
40        else compressionLevel= file->compression_level.isEmpty() ? 0 :file->compression_level ;
41      }
42
43      CNc4DataOutput::~CNc4DataOutput(void)
44    { /* Ne rien faire de plus */ }
45
46      ///--------------------------------------------------------------
47
48      const StdString & CNc4DataOutput::getFileName(void) const
49      {
50         return (this->filename);
51      }
52
53      //---------------------------------------------------------------
54
55      void CNc4DataOutput::writeDomain_(CDomain* domain)
56      TRY
57      {
58        StdString lonName,latName ;
59
60        if (domain->type == CDomain::type_attr::unstructured)
61        {
62          if (SuperClassWriter::useCFConvention)
63            writeUnstructuredDomain(domain) ;
64          else
65            writeUnstructuredDomainUgrid(domain) ;
66          return ;
67        }
68
69         CContext* context = CContext::getCurrent() ;
70         if (domain->IsWritten(this->filename)) return;
71         domain->checkAttributes();
72
73         if (domain->isEmpty())
74           if (SuperClass::type==MULTI_FILE) return;
75
76
77         // Check that the name associated to the current element is not in conflict with an existing element (due to CGrid::duplicateSentGrid)
78         if (!domain->lonvalue.isEmpty() )
79         {
80           // The hash of the element will be associated to the default element name (= map key), and to the name really written
81           int globalHash = domain->computeAttributesHash( comm_file ); // Need a MPI_Comm to distribute without redundancy some attributs (value)
82
83           StdString defaultNameKey = domain->getDomainOutputName();
84           if ( !relDomains_.count ( defaultNameKey ) )
85           {
86             // if defaultNameKey not in the map, write the element such as it is defined
87             relDomains_.insert( make_pair( defaultNameKey, make_pair(globalHash, domain) ) );
88           }
89           else // look if a hash associated this key is equal
90           {
91             bool elementIsInMap(false);
92             auto defaultNameKeyElements = relDomains_.equal_range( defaultNameKey );
93             for (auto it = defaultNameKeyElements.first; it != defaultNameKeyElements.second; it++)
94             {
95               if ( it->second.first == globalHash )
96               {
97                 // if yes, associate the same ids to current element
98                 domain->name = it->second.second->getDomainOutputName();
99                 domain->lon_name = it->second.second->lon_name;
100                 domain->lat_name = it->second.second->lat_name;
101                 elementIsInMap = true;
102               }
103             }
104             // if no : inheritance has been excessive, define new names and store it (could be used by another grid)
105             if (!elementIsInMap)  // ! in MAP
106             {
107               domain->name =  domain->getId();
108               domain->lon_name = "lon_"+domain->getId();
109               domain->lat_name = "lat_"+domain->getId();
110               relDomains_.insert( make_pair( defaultNameKey, make_pair(globalHash, domain) ) ) ;         
111             }
112           }
113         }
114
115         
116         std::vector<StdString> dim0, dim1;
117         StdString domid = domain->getDomainOutputName();
118         StdString appendDomid  = (singleDomain) ? "" : "_"+domid ;
119         if (isWrittenDomain(domid)) return ;
120         else setWrittenDomain(domid);
121       
122         int nvertex = (domain->nvertex.isEmpty()) ? 0 : domain->nvertex;
123
124
125        StdString dimXid, dimYid ;
126
127        nc_type typePrec ;
128        if (domain->prec.isEmpty()) typePrec =  NC_FLOAT ;
129        else if (domain->prec==4)  typePrec =  NC_FLOAT ;
130        else if (domain->prec==8)   typePrec =  NC_DOUBLE ;
131         
132         bool isRegularDomain = (domain->type == CDomain::type_attr::rectilinear);
133         switch (domain->type)
134         {
135           case CDomain::type_attr::curvilinear :
136
137             if (domain->lon_name.isEmpty()) lonName = "nav_lon";
138             else lonName = domain->lon_name;
139
140             if (domain->lat_name.isEmpty()) latName = "nav_lat";
141             else latName = domain->lat_name;
142
143             if (domain->dim_i_name.isEmpty()) dimXid=StdString("x").append(appendDomid);
144             else dimXid=domain->dim_i_name.getValue() + appendDomid;
145
146             if (domain->dim_j_name.isEmpty()) dimYid=StdString("y").append(appendDomid);
147             else dimYid=domain->dim_j_name.getValue() + appendDomid;
148
149             break ;
150
151           case CDomain::type_attr::rectilinear :
152
153             if (domain->lon_name.isEmpty())
154             {
155               if (domain->dim_i_name.isEmpty())
156                   lonName = "lon";
157               else
158                 lonName = domain->dim_i_name.getValue();
159             }
160             else lonName = domain->lon_name;
161
162             if (domain->lat_name.isEmpty())
163             {
164               if (domain->dim_j_name.isEmpty())
165                 latName = "lat";
166               else
167                 latName = domain->dim_j_name.getValue();
168             }
169             else latName = domain->lat_name;
170             
171             if (domain->dim_i_name.isEmpty()) dimXid = lonName+appendDomid;
172             else dimXid = domain->dim_i_name.getValue()+appendDomid;
173
174             if (domain->dim_j_name.isEmpty()) dimYid = latName+appendDomid;
175             else dimYid = domain->dim_j_name.getValue()+appendDomid;
176             break;
177         }
178
179         StdString dimVertId = StdString("nvertex").append(appendDomid);
180
181         string lonid,latid,bounds_lonid,bounds_latid ;
182         string areaId = "area" + appendDomid;
183
184         try
185         {
186           switch (SuperClass::type)
187           {
188              case (MULTI_FILE) :
189              {
190                 switch (domain->type)
191                 {
192                   case CDomain::type_attr::curvilinear :
193                     dim0.push_back(dimYid); dim0.push_back(dimXid);
194                     lonid = lonName+appendDomid;
195                     latid = latName+appendDomid;
196                     break ;
197                   case CDomain::type_attr::rectilinear :
198                     lonid = lonName+appendDomid;
199                     latid = latName+appendDomid;
200                     dim0.push_back(dimYid);
201                     dim1.push_back(dimXid);
202                     break;
203                 }
204                 if (!domain->bounds_lon_name.isEmpty()) bounds_lonid = domain->bounds_lon_name;
205                 else bounds_lonid = "bounds_"+lonName+appendDomid;
206                 if (!domain->bounds_lat_name.isEmpty()) bounds_latid = domain->bounds_lat_name;
207                 else bounds_latid = "bounds_"+latName+appendDomid;
208
209                 SuperClassWriter::addDimension(dimXid, domain->ni);
210                 SuperClassWriter::addDimension(dimYid, domain->nj);
211
212                 if (domain->hasBounds)
213                   SuperClassWriter::addDimension(dimVertId, domain->nvertex);
214
215                 if (context->intraCommSize_ > 1)
216                 {
217                   this->writeLocalAttributes(domain->ibegin,
218                                              domain->ni,
219                                              domain->jbegin,
220                                              domain->nj,
221                                              appendDomid);
222
223                   if (singleDomain)
224                    this->writeLocalAttributes_IOIPSL(dimXid, dimYid,
225                                                      domain->ibegin,
226                                                      domain->ni,
227                                                      domain->jbegin,
228                                                      domain->nj,
229                                                      domain->ni_glo,domain->nj_glo,
230                                                      context->intraCommRank_,context->intraCommSize_);
231                 }
232
233                 if (domain->hasLonLat)
234                 {
235                   switch (domain->type)
236                   {
237                     case CDomain::type_attr::curvilinear :
238                       SuperClassWriter::addVariable(latid, typePrec, dim0, compressionLevel);
239                       SuperClassWriter::addVariable(lonid, typePrec, dim0, compressionLevel);
240                       break ;
241                      case CDomain::type_attr::rectilinear :
242                        SuperClassWriter::addVariable(latid, typePrec, dim0, compressionLevel);
243                        SuperClassWriter::addVariable(lonid, typePrec, dim1, compressionLevel);
244                        break ;
245                   }
246
247                   this->writeAxisAttributes(lonid, isRegularDomain ? "X" : "", "longitude", "Longitude", "degrees_east", domid);
248                   this->writeAxisAttributes(latid, isRegularDomain ? "Y" : "", "latitude", "Latitude", "degrees_north", domid);
249
250                   if (domain->hasBounds)
251                   {
252                     SuperClassWriter::addAttribute("bounds", bounds_lonid, &lonid);
253                     SuperClassWriter::addAttribute("bounds", bounds_latid, &latid);
254
255                     dim0.clear();
256                     dim0.push_back(dimYid);
257                     dim0.push_back(dimXid);
258                     dim0.push_back(dimVertId);
259                     SuperClassWriter::addVariable(bounds_lonid, typePrec, dim0, compressionLevel);
260                     SuperClassWriter::addVariable(bounds_latid, typePrec, dim0, compressionLevel);
261                   }
262                 }
263
264                 dim0.clear();
265                 dim0.push_back(dimYid);
266                 dim0.push_back(dimXid);
267
268                 if (domain->hasArea)
269                 {
270                   SuperClassWriter::addVariable(areaId, typePrec, dim0, compressionLevel);
271                   SuperClassWriter::addAttribute("standard_name", StdString("cell_area"), &areaId);
272                   SuperClassWriter::addAttribute("units", StdString("m2"), &areaId);
273                 }
274
275                 SuperClassWriter::definition_end();
276
277                 if (domain->hasLonLat)
278                 {
279                   switch (domain->type)
280                   {
281                     case CDomain::type_attr::curvilinear :                       
282                       SuperClassWriter::writeData(domain->latvalue, latid, isCollective, 0);
283                       SuperClassWriter::writeData(domain->lonvalue, lonid, isCollective, 0);
284                       break;
285                     case CDomain::type_attr::rectilinear :
286                       CArray<double,1> lat = domain->latvalue(Range(fromStart,toEnd,domain->ni)) ;
287                       SuperClassWriter::writeData(CArray<double,1>(lat.copy()), latid, isCollective, 0);
288                       CArray<double,1> lon = domain->lonvalue(Range(0,domain->ni-1)) ;
289                       SuperClassWriter::writeData(CArray<double,1>(lon.copy()), lonid, isCollective, 0);
290                       break;
291                   }
292
293                   if (domain->hasBounds)
294                   {
295                     SuperClassWriter::writeData(domain->bounds_lonvalue, bounds_lonid, isCollective, 0);
296                     SuperClassWriter::writeData(domain->bounds_latvalue, bounds_latid, isCollective, 0);
297                   }
298                 }
299
300                 if (domain->hasArea)
301                 {
302                   SuperClassWriter::writeData(domain->areavalue, areaId, isCollective, 0);                   
303                 }
304
305                 SuperClassWriter::definition_start();
306
307                 break;
308              }
309              case (ONE_FILE) :
310              {
311                SuperClassWriter::addDimension(dimXid, domain->ni_glo);
312                SuperClassWriter::addDimension(dimYid, domain->nj_glo);
313
314                 if (domain->hasBounds)
315                   SuperClassWriter::addDimension(dimVertId, domain->nvertex);
316
317                 if (domain->hasLonLat)
318                 {
319                   switch (domain->type)
320                   {
321                     case CDomain::type_attr::curvilinear :
322                       dim0.push_back(dimYid); dim0.push_back(dimXid);
323                       lonid = lonName+appendDomid;
324                       latid = latName+appendDomid;
325                       SuperClassWriter::addVariable(latid, typePrec, dim0, compressionLevel);
326                       SuperClassWriter::addVariable(lonid, typePrec, dim0, compressionLevel);
327                       break;
328
329                     case CDomain::type_attr::rectilinear :
330                       dim0.push_back(dimYid);
331                       dim1.push_back(dimXid);
332                       lonid = lonName+appendDomid;
333                       latid = latName+appendDomid;
334                       SuperClassWriter::addVariable(latid, typePrec, dim0, compressionLevel);
335                       SuperClassWriter::addVariable(lonid, typePrec, dim1, compressionLevel);
336                       break;
337                   }
338                   if (!domain->bounds_lon_name.isEmpty()) bounds_lonid = domain->bounds_lon_name;
339                   else bounds_lonid = "bounds_"+lonName+appendDomid;
340                   if (!domain->bounds_lat_name.isEmpty()) bounds_latid = domain->bounds_lat_name;
341                   else bounds_latid = "bounds_"+latName+appendDomid;
342
343                   this->writeAxisAttributes
344                      (lonid, isRegularDomain ? "X" : "", "longitude", "Longitude", "degrees_east", domid);
345                   this->writeAxisAttributes
346                      (latid, isRegularDomain ? "Y" : "", "latitude", "Latitude", "degrees_north", domid);
347
348                   if (domain->hasBounds)
349                   {
350                     SuperClassWriter::addAttribute("bounds", bounds_lonid, &lonid);
351                     SuperClassWriter::addAttribute("bounds", bounds_latid, &latid);
352
353                     dim0.clear();
354                     dim0.push_back(dimYid);
355                     dim0.push_back(dimXid);
356                     dim0.push_back(dimVertId);
357                     SuperClassWriter::addVariable(bounds_lonid, typePrec, dim0, compressionLevel);
358                     SuperClassWriter::addVariable(bounds_latid, typePrec, dim0, compressionLevel);
359                   }
360                 }
361
362                 if (domain->hasArea)
363                 {
364                   dim0.clear();
365                   dim0.push_back(dimYid); dim0.push_back(dimXid);
366                   SuperClassWriter::addVariable(areaId, typePrec, dim0, compressionLevel);
367                   SuperClassWriter::addAttribute("standard_name", StdString("cell_area"), &areaId);
368                   SuperClassWriter::addAttribute("units", StdString("m2"), &areaId);
369                   dim0.clear();
370                 }
371
372                 SuperClassWriter::definition_end();
373
374                 switch (domain->type)
375                 {
376                   case CDomain::type_attr::curvilinear :
377                   {
378                     std::vector<StdSize> start(2) ;
379                     std::vector<StdSize> count(2) ;
380                     start[1]=domain->ibegin;
381                     start[0]=domain->jbegin;
382                     count[1]=domain->ni ; count[0]=domain->nj ;
383
384                     if (domain->hasLonLat)
385                     {
386                       SuperClassWriter::writeData(domain->latvalue, latid, isCollective, 0,&start,&count);
387                       SuperClassWriter::writeData(domain->lonvalue, lonid, isCollective, 0,&start,&count);
388                     }
389                     break;
390                   }
391                   case CDomain::type_attr::rectilinear :
392                   {
393                     if (domain->hasLonLat)
394                     {
395                       std::vector<StdSize> start(1) ;
396                       std::vector<StdSize> count(1) ;
397                       
398                       start[0]=domain->jbegin;
399                       count[0]=domain->nj;
400                       CArray<double,1> lat = domain->latvalue(Range(fromStart,toEnd,domain->ni));
401                       SuperClassWriter::writeData(CArray<double,1>(lat.copy()), latid, isCollective, 0,&start,&count);
402
403                       start[0]=domain->ibegin;
404                       count[0]=domain->ni;
405                       CArray<double,1> lon = domain->lonvalue(Range(0,domain->ni-1));
406                       SuperClassWriter::writeData(CArray<double,1>(lon.copy()), lonid, isCollective, 0,&start,&count);
407                     }
408                     break;
409                   }
410                 }
411
412                 if (domain->hasBounds)
413                 {
414                   std::vector<StdSize> start(3);
415                   std::vector<StdSize> count(3);
416                   if (domain->isEmpty())
417                   {
418                     start[2] = start[1] = start[0] = 0;
419                     count[2] = count[1] = count[0] = 0;
420                   }
421                   else
422                   {
423                     start[2] = 0;
424                     start[1] = domain->ibegin;
425                     start[0] = domain->jbegin;
426                     count[2] = domain->nvertex;
427                     count[1] = domain->ni;
428                     count[0] = domain->nj;
429                   }
430                 
431                   SuperClassWriter::writeData(domain->bounds_lonvalue, bounds_lonid, isCollective, 0, &start, &count); // will probably not working for rectilinear
432                   SuperClassWriter::writeData(domain->bounds_latvalue, bounds_latid, isCollective, 0, &start, &count);
433                 }
434
435                 if (domain->hasArea)
436                 {
437                   std::vector<StdSize> start(2);
438                   std::vector<StdSize> count(2);
439
440                   start[1] = domain->ibegin;
441                   start[0] = domain->jbegin;
442                   count[1] = domain->ni;
443                   count[0] = domain->nj;
444                   
445                   SuperClassWriter::writeData(domain->areavalue, areaId, isCollective, 0, &start, &count);
446                 }
447
448                 SuperClassWriter::definition_start();
449                 break;
450              }
451              default :
452                 ERROR("CNc4DataOutput::writeDomain(domain)",
453                       << "[ type = " << SuperClass::type << "]"
454                       << " not implemented yet !");
455           }
456         }
457         catch (CNetCdfException& e)
458         {
459           StdString msg("On writing the domain : ");
460           msg.append(domid); msg.append("\n");
461           msg.append("In the context : ");
462           msg.append(context->getId()); msg.append("\n");
463           msg.append(e.what());
464           ERROR("CNc4DataOutput::writeDomain_(CDomain* domain)", << msg);
465         }
466
467         domain->addRelFile(this->filename);
468      }
469      CATCH
470
471    //--------------------------------------------------------------
472
473    void CNc4DataOutput::writeUnstructuredDomainUgrid(CDomain* domain)
474    {
475      CContext* context = CContext::getCurrent() ;
476
477      if (domain->IsWritten(this->filename)) return;
478
479      StdString domid = domain->getDomainOutputName();
480
481      // The first domain for the same mesh that will be written is that with the highest value of nvertex.
482      // Thus the entire mesh connectivity will be generated at once.
483      if (isWrittenDomain(domid)) return ;
484      else setWrittenDomain(domid);
485
486      domain->checkAttributes();
487      if (domain->isEmpty())
488        if (SuperClass::type==MULTI_FILE) return ;
489
490     nc_type typePrec ;
491     if (domain->prec.isEmpty()) typePrec =  NC_FLOAT ;
492     else if (domain->prec==4)  typePrec =  NC_FLOAT ;
493     else if (domain->prec==8)   typePrec =  NC_DOUBLE ;
494
495      std::vector<StdString> dim0;
496      StdString domainName = domain->name;
497      domain->assignMesh(domainName, domain->nvertex);
498      domain->mesh->createMeshEpsilon(context->intraComm_, domain->lonvalue, domain->latvalue, domain->bounds_lonvalue, domain->bounds_latvalue);
499
500      StdString node_x = domainName + "_node_x";
501      StdString node_y = domainName + "_node_y";
502
503      StdString edge_x = domainName + "_edge_x";
504      StdString edge_y = domainName + "_edge_y";
505      StdString edge_nodes = domainName + "_edge_nodes";
506
507      StdString face_x = domainName + "_face_x";
508      StdString face_y = domainName + "_face_y";
509      StdString face_nodes = domainName + "_face_nodes";
510      StdString face_edges = domainName + "_face_edges";
511      StdString edge_faces = domainName + "_edge_face_links";
512      StdString face_faces = domainName + "_face_links";
513
514      StdString dimNode = "n" + domainName + "_node";
515      StdString dimEdge = "n" + domainName + "_edge";
516      StdString dimFace = "n" + domainName + "_face";
517      StdString dimVertex = "n" + domainName + "_vertex";
518      StdString dimTwo = "Two";
519
520      if (!SuperClassWriter::dimExist(dimTwo)) SuperClassWriter::addDimension(dimTwo, 2);
521      dim0.clear();
522      SuperClassWriter::addVariable(domainName, NC_INT, dim0, compressionLevel);
523      SuperClassWriter::addAttribute("cf_role", StdString("mesh_topology"), &domainName);
524      SuperClassWriter::addAttribute("long_name", StdString("Topology data of 2D unstructured mesh"), &domainName);
525      SuperClassWriter::addAttribute("topology_dimension", 2, &domainName);
526      SuperClassWriter::addAttribute("node_coordinates", node_x + " " + node_y, &domainName);
527
528      try
529      {
530        switch (SuperClass::type)
531        {
532          case (ONE_FILE) :
533          {
534            // Adding nodes
535            if (domain->nvertex == 1)
536            {
537              if (!SuperClassWriter::varExist(node_x) || !SuperClassWriter::varExist(node_y))
538              {
539                SuperClassWriter::addDimension(dimNode, domain->ni_glo);
540                dim0.clear();
541                dim0.push_back(dimNode);
542                SuperClassWriter::addVariable(node_x, typePrec, dim0, compressionLevel);
543                SuperClassWriter::addAttribute("standard_name", StdString("longitude"), &node_x);
544                SuperClassWriter::addAttribute("long_name", StdString("Longitude of mesh nodes."), &node_x);
545                SuperClassWriter::addAttribute("units", StdString("degrees_east"), &node_x);
546                SuperClassWriter::addVariable(node_y, typePrec, dim0, compressionLevel);
547                SuperClassWriter::addAttribute("standard_name", StdString("latitude"), &node_y);
548                SuperClassWriter::addAttribute("long_name", StdString("Latitude of mesh nodes."), &node_y);
549                SuperClassWriter::addAttribute("units", StdString("degrees_north"), &node_y);
550              }
551            } // domain->nvertex == 1
552
553            // Adding edges and nodes, if nodes have not been defined previously
554            if (domain->nvertex == 2)
555            {
556              if (!SuperClassWriter::varExist(node_x) || !SuperClassWriter::varExist(node_y))
557              {
558                SuperClassWriter::addDimension(dimNode, domain->mesh->nbNodesGlo);
559                dim0.clear();
560                dim0.push_back(dimNode);
561                SuperClassWriter::addVariable(node_x, typePrec, dim0, compressionLevel);
562                SuperClassWriter::addAttribute("standard_name", StdString("longitude"), &node_x);
563                SuperClassWriter::addAttribute("long_name", StdString("Longitude of mesh nodes."), &node_x);
564                SuperClassWriter::addAttribute("units", StdString("degrees_east"), &node_x);
565                SuperClassWriter::addVariable(node_y, typePrec, dim0, compressionLevel);
566                SuperClassWriter::addAttribute("standard_name", StdString("latitude"), &node_y);
567                SuperClassWriter::addAttribute("long_name", StdString("Latitude of mesh nodes."), &node_y);
568                SuperClassWriter::addAttribute("units", StdString("degrees_north"), &node_y);
569              }
570              SuperClassWriter::addAttribute("edge_node_connectivity", edge_nodes, &domainName);
571              SuperClassWriter::addAttribute("edge_coordinates", edge_x + " " + edge_y, &domainName);
572              SuperClassWriter::addDimension(dimEdge, domain->ni_glo);
573              dim0.clear();
574              dim0.push_back(dimEdge);
575              SuperClassWriter::addVariable(edge_x, typePrec, dim0, compressionLevel);
576              SuperClassWriter::addAttribute("standard_name", StdString("longitude"), &edge_x);
577              SuperClassWriter::addAttribute("long_name", StdString("Characteristic longitude of mesh edges."), &edge_x);
578              SuperClassWriter::addAttribute("units", StdString("degrees_east"), &edge_x);
579              SuperClassWriter::addVariable(edge_y, typePrec, dim0, compressionLevel);
580              SuperClassWriter::addAttribute("standard_name", StdString("latitude"), &edge_y);
581              SuperClassWriter::addAttribute("long_name", StdString("Characteristic latitude of mesh edges."), &edge_y);
582              SuperClassWriter::addAttribute("units", StdString("degrees_north"), &edge_y);
583              dim0.clear();
584              dim0.push_back(dimEdge);
585              dim0.push_back(dimTwo);
586              SuperClassWriter::addVariable(edge_nodes, NC_INT, dim0, compressionLevel);
587              SuperClassWriter::addAttribute("cf_role", StdString("edge_node_connectivity"), &edge_nodes);
588              SuperClassWriter::addAttribute("long_name", StdString("Maps every edge/link to two nodes that it connects."), &edge_nodes);
589              SuperClassWriter::addAttribute("start_index", 0, &edge_nodes);
590            } // domain->nvertex == 2
591
592            // Adding faces, edges, and nodes, if edges and nodes have not been defined previously
593            if (domain->nvertex > 2)
594            {
595              // Nodes
596              if (!SuperClassWriter::varExist(node_x) || !SuperClassWriter::varExist(node_y))
597              {
598                SuperClassWriter::addDimension(dimNode, domain->mesh->nbNodesGlo);
599                dim0.clear();
600                dim0.push_back(dimNode);
601                SuperClassWriter::addVariable(node_x, typePrec, dim0, compressionLevel);
602                SuperClassWriter::addAttribute("standard_name", StdString("longitude"), &node_x);
603                SuperClassWriter::addAttribute("long_name", StdString("Longitude of mesh nodes."), &node_x);
604                SuperClassWriter::addAttribute("units", StdString("degrees_east"), &node_x);
605                SuperClassWriter::addVariable(node_y, typePrec, dim0, compressionLevel);
606                SuperClassWriter::addAttribute("standard_name", StdString("latitude"), &node_y);
607                SuperClassWriter::addAttribute("long_name", StdString("Latitude of mesh nodes."), &node_y);
608                SuperClassWriter::addAttribute("units", StdString("degrees_north"), &node_y);
609              }
610              if (!SuperClassWriter::varExist(edge_x) || !SuperClassWriter::varExist(edge_y))
611              {
612                SuperClassWriter::addAttribute("edge_coordinates", edge_x + " " + edge_y, &domainName);
613                SuperClassWriter::addAttribute("edge_node_connectivity", edge_nodes, &domainName);
614                SuperClassWriter::addDimension(dimEdge, domain->mesh->nbEdgesGlo);
615                dim0.clear();
616                dim0.push_back(dimEdge);
617                SuperClassWriter::addVariable(edge_x, typePrec, dim0, compressionLevel);
618                SuperClassWriter::addAttribute("standard_name", StdString("longitude"), &edge_x);
619                SuperClassWriter::addAttribute("long_name", StdString("Characteristic longitude of mesh edges."), &edge_x);
620                SuperClassWriter::addAttribute("units", StdString("degrees_east"), &edge_x);
621                SuperClassWriter::addVariable(edge_y, typePrec, dim0, compressionLevel);
622                SuperClassWriter::addAttribute("standard_name", StdString("latitude"), &edge_y);
623                SuperClassWriter::addAttribute("long_name", StdString("Characteristic latitude of mesh edges."), &edge_y);
624                SuperClassWriter::addAttribute("units", StdString("degrees_north"), &edge_y);
625                dim0.clear();
626                dim0.push_back(dimEdge);
627                dim0.push_back(dimTwo);
628                SuperClassWriter::addVariable(edge_nodes, NC_INT, dim0, compressionLevel);
629                SuperClassWriter::addAttribute("cf_role", StdString("edge_node_connectivity"), &edge_nodes);
630                SuperClassWriter::addAttribute("long_name", StdString("Maps every edge/link to two nodes that it connects."), &edge_nodes);
631                SuperClassWriter::addAttribute("start_index", 0, &edge_nodes);
632              }
633              SuperClassWriter::addAttribute("face_coordinates", face_x + " " + face_y, &domainName);
634              SuperClassWriter::addAttribute("face_node_connectivity", face_nodes, &domainName);
635              SuperClassWriter::addDimension(dimFace, domain->ni_glo);
636              SuperClassWriter::addDimension(dimVertex, domain->nvertex);
637              dim0.clear();
638              dim0.push_back(dimFace);
639              SuperClassWriter::addVariable(face_x, typePrec, dim0, compressionLevel);
640              SuperClassWriter::addAttribute("standard_name", StdString("longitude"), &face_x);
641              SuperClassWriter::addAttribute("long_name", StdString("Characteristic longitude of mesh faces."), &face_x);
642              SuperClassWriter::addAttribute("units", StdString("degrees_east"), &face_x);
643              SuperClassWriter::addVariable(face_y, typePrec, dim0, compressionLevel);
644              SuperClassWriter::addAttribute("standard_name", StdString("latitude"), &face_y);
645              SuperClassWriter::addAttribute("long_name", StdString("Characteristic latitude of mesh faces."), &face_y);
646              SuperClassWriter::addAttribute("units", StdString("degrees_north"), &face_y);
647              dim0.clear();
648              dim0.push_back(dimFace);
649              dim0.push_back(dimVertex);
650              SuperClassWriter::addVariable(face_nodes, NC_INT, dim0, compressionLevel);
651              SuperClassWriter::addAttribute("cf_role", StdString("face_node_connectivity"), &face_nodes);
652              SuperClassWriter::addAttribute("long_name", StdString("Maps every face to its corner nodes."), &face_nodes);
653              SuperClassWriter::addAttribute("start_index", 0, &face_nodes);
654              dim0.clear();
655              dim0.push_back(dimFace);
656              dim0.push_back(dimVertex);
657              SuperClassWriter::addVariable(face_edges, NC_INT, dim0, compressionLevel);
658              SuperClassWriter::addAttribute("cf_role", StdString("face_edge_connectivity"), &face_edges);
659              SuperClassWriter::addAttribute("long_name", StdString("Maps every face to its edges."), &face_edges);
660              SuperClassWriter::addAttribute("start_index", 0, &face_edges);
661              SuperClassWriter::addAttribute("_FillValue", 999999, &face_edges);
662              dim0.clear();
663              dim0.push_back(dimEdge);
664              dim0.push_back(dimTwo);
665              SuperClassWriter::addVariable(edge_faces, NC_INT, dim0, compressionLevel);
666              SuperClassWriter::addAttribute("cf_role", StdString("edge_face connectivity"), &edge_faces);
667              SuperClassWriter::addAttribute("long_name", StdString("neighbor faces for edges"), &edge_faces);
668              SuperClassWriter::addAttribute("start_index", 0, &edge_faces);
669              SuperClassWriter::addAttribute("_FillValue", -999, &edge_faces);
670              SuperClassWriter::addAttribute("comment", StdString("missing neighbor faces are indicated using _FillValue"), &edge_faces);
671              dim0.clear();
672              dim0.push_back(dimFace);
673              dim0.push_back(dimVertex);
674              SuperClassWriter::addVariable(face_faces, NC_INT, dim0, compressionLevel);
675              SuperClassWriter::addAttribute("cf_role", StdString("face_face connectivity"), &face_faces);
676              SuperClassWriter::addAttribute("long_name", StdString("Indicates which other faces neighbor each face"), &face_faces);
677              SuperClassWriter::addAttribute("start_index", 0, &face_faces);
678              SuperClassWriter::addAttribute("_FillValue", 999999, &face_faces);
679              SuperClassWriter::addAttribute("flag_values", -1, &face_faces);
680              SuperClassWriter::addAttribute("flag_meanings", StdString("out_of_mesh"), &face_faces);
681            } // domain->nvertex > 2
682
683            SuperClassWriter::definition_end();
684
685            std::vector<StdSize> startEdges(1) ;
686            std::vector<StdSize> countEdges(1) ;
687            std::vector<StdSize> startNodes(1) ;
688            std::vector<StdSize> countNodes(1) ;
689            std::vector<StdSize> startFaces(1) ;
690            std::vector<StdSize> countFaces(1) ;
691            std::vector<StdSize> startEdgeNodes(2) ;
692            std::vector<StdSize> countEdgeNodes(2) ;
693            std::vector<StdSize> startEdgeFaces(2) ;
694            std::vector<StdSize> countEdgeFaces(2) ;
695            std::vector<StdSize> startFaceConctv(2) ;
696            std::vector<StdSize> countFaceConctv(2) ;
697
698            if (domain->nvertex == 1)
699            {
700              if (domain->isEmpty())
701               {
702                 startNodes[0]=0 ;
703                 countNodes[0]=0 ;
704               }
705               else
706               {
707                 startNodes[0] = domain->ibegin;
708                 countNodes[0] = domain->ni ;
709               }
710
711              SuperClassWriter::writeData(domain->mesh->node_lat, node_y, isCollective, 0, &startNodes, &countNodes);
712              SuperClassWriter::writeData(domain->mesh->node_lon, node_x, isCollective, 0, &startNodes, &countNodes);
713            }
714            else if (domain->nvertex == 2)
715            {
716              if (domain->isEmpty())
717               {
718                startEdges[0]=0 ;
719                countEdges[0]=0 ;
720                startNodes[0]=0 ;
721                countNodes[0]=0 ;
722                startEdgeNodes[0]=0;
723                startEdgeNodes[1]=0;
724                countEdgeNodes[0]=0;
725                countEdgeNodes[1]=0;
726
727               }
728               else
729               {
730                 startEdges[0] = domain->ibegin;
731                 countEdges[0] = domain->ni;
732                 startNodes[0] = domain->mesh->node_start;
733                 countNodes[0] = domain->mesh->node_count;
734                 startEdgeNodes[0] = domain->ibegin;
735                 startEdgeNodes[1] = 0;
736                 countEdgeNodes[0] = domain->ni;
737                 countEdgeNodes[1] = 2;
738               }
739              SuperClassWriter::writeData(domain->mesh->node_lat, node_y, isCollective, 0, &startNodes, &countNodes);
740              SuperClassWriter::writeData(domain->mesh->node_lon, node_x, isCollective, 0, &startNodes, &countNodes);
741              SuperClassWriter::writeData(domain->mesh->edge_lat, edge_y, isCollective, 0, &startEdges, &countEdges);
742              SuperClassWriter::writeData(domain->mesh->edge_lon, edge_x, isCollective, 0, &startEdges, &countEdges);
743              SuperClassWriter::writeData(domain->mesh->edge_nodes, edge_nodes, isCollective, 0, &startEdgeNodes, &countEdgeNodes);
744            }
745            else
746            {
747              if (domain->isEmpty())
748               {
749                 startFaces[0] = 0 ;
750                 countFaces[0] = 0 ;
751                 startNodes[0] = 0;
752                 countNodes[0] = 0;
753                 startEdges[0] = 0;
754                 countEdges[0] = 0;
755                 startEdgeFaces[0] = 0;
756                 startEdgeFaces[1] = 0;
757                 countEdgeFaces[0] = 0;
758                 countEdgeFaces[1] = 0;
759                 startFaceConctv[0] = 0;
760                 startFaceConctv[1] = 0;
761                 countFaceConctv[0] = 0;
762                 countFaceConctv[1] = 0;
763               }
764               else
765               {
766                 startFaces[0] = domain->ibegin;
767                 countFaces[0] = domain->ni ;
768                 startNodes[0] = domain->mesh->node_start;
769                 countNodes[0] = domain->mesh->node_count;
770                 startEdges[0] = domain->mesh->edge_start;
771                 countEdges[0] = domain->mesh->edge_count;
772                 startEdgeNodes[0] = domain->mesh->edge_start;
773                 startEdgeNodes[1] = 0;
774                 countEdgeNodes[0] = domain->mesh->edge_count;
775                 countEdgeNodes[1]= 2;
776                 startEdgeFaces[0] = domain->mesh->edge_start;
777                 startEdgeFaces[1]= 0;
778                 countEdgeFaces[0] = domain->mesh->edge_count;
779                 countEdgeFaces[1]= 2;
780                 startFaceConctv[0] = domain->ibegin;
781                 startFaceConctv[1] = 0;
782                 countFaceConctv[0] = domain->ni;
783                 countFaceConctv[1] = domain->nvertex;
784               }
785              SuperClassWriter::writeData(domain->mesh->node_lat, node_y, isCollective, 0, &startNodes, &countNodes);
786              SuperClassWriter::writeData(domain->mesh->node_lon, node_x, isCollective, 0, &startNodes, &countNodes);
787              SuperClassWriter::writeData(domain->mesh->edge_lat, edge_y, isCollective, 0, &startEdges, &countEdges);
788              SuperClassWriter::writeData(domain->mesh->edge_lon, edge_x, isCollective, 0, &startEdges, &countEdges);
789              SuperClassWriter::writeData(domain->mesh->edge_nodes, edge_nodes, isCollective, 0, &startEdgeNodes, &countEdgeNodes);
790              SuperClassWriter::writeData(domain->mesh->face_lat, face_y, isCollective, 0, &startFaces, &countFaces);
791              SuperClassWriter::writeData(domain->mesh->face_lon, face_x, isCollective, 0, &startFaces, &countFaces);
792              SuperClassWriter::writeData(domain->mesh->face_nodes, face_nodes, isCollective, 0, &startFaceConctv, &countFaceConctv);
793              SuperClassWriter::writeData(domain->mesh->face_edges, face_edges, isCollective, 0, &startFaceConctv, &countFaceConctv);
794              SuperClassWriter::writeData(domain->mesh->edge_faces, edge_faces, isCollective, 0, &startEdgeFaces, &countEdgeFaces);
795              SuperClassWriter::writeData(domain->mesh->face_faces, face_faces, isCollective, 0, &startFaceConctv, &countFaceConctv);
796            }
797            SuperClassWriter::definition_start();
798
799            break;
800          } // ONE_FILE
801
802          case (MULTI_FILE) :
803          {
804            ERROR("CNc4DataOutput::writeDomain(domain)",
805            << "[ type = multiple_file ]"
806            << " is not yet implemented for UGRID files !");
807            break;
808          }
809
810          default :
811          ERROR("CNc4DataOutput::writeDomain(domain)",
812          << "[ type = " << SuperClass::type << "]"
813          << " not implemented yet !");
814          } // switch
815        } // try
816
817        catch (CNetCdfException& e)
818        {
819          StdString msg("On writing the domain : ");
820          msg.append(domid); msg.append("\n");
821          msg.append("In the context : ");
822          msg.append(context->getId()); msg.append("\n");
823          msg.append(e.what());
824          ERROR("CNc4DataOutput::writeUnstructuredDomainUgrid(CDomain* domain)", << msg);
825        }
826
827  domain->addRelFile(this->filename);
828  }
829
830    //--------------------------------------------------------------
831
832    void CNc4DataOutput::writeUnstructuredDomain(CDomain* domain)
833      {
834         CContext* context = CContext::getCurrent() ;
835
836         if (domain->IsWritten(this->filename)) return;
837         domain->checkAttributes();
838
839         if (domain->isEmpty())
840           if (SuperClass::type==MULTI_FILE) return ;
841
842         std::vector<StdString> dim0, dim1;
843         StdString domid = domain->getDomainOutputName();
844         if (isWrittenDomain(domid)) return ;
845         else setWrittenDomain(domid);
846
847         StdString appendDomid  = (singleDomain) ? "" : "_"+domid ;
848
849         StdString lonName,latName, cellName ;
850         if (domain->lon_name.isEmpty()) lonName = "lon";
851         else lonName = domain->lon_name;
852
853         if (domain->lat_name.isEmpty()) latName = "lat";
854         else latName = domain->lat_name;
855
856         if (!domain->dim_i_name.isEmpty()) cellName=domain->dim_i_name;
857         else cellName="cell";
858         StdString dimXid = cellName+appendDomid;
859         StdString dimVertId = StdString("nvertex").append(appendDomid);
860
861         string lonid,latid,bounds_lonid,bounds_latid ;
862         string areaId = "area" + appendDomid;
863
864         nc_type typePrec ;
865         if (domain->prec.isEmpty()) typePrec =  NC_FLOAT ;
866         else if (domain->prec==4)  typePrec =  NC_FLOAT ;
867         else if (domain->prec==8)   typePrec =  NC_DOUBLE ;
868
869         int nvertex = (domain->nvertex.isEmpty()) ? 0 : domain->nvertex;
870
871         try
872         {
873           switch (SuperClass::type)
874           {
875              case (MULTI_FILE) :
876              {
877                 dim0.push_back(dimXid);
878                 SuperClassWriter::addDimension(dimXid, domain->ni);
879
880                 lonid = lonName+appendDomid;
881                 latid = latName+appendDomid;
882                 if (!domain->bounds_lon_name.isEmpty()) bounds_lonid = domain->bounds_lon_name;
883                 else bounds_lonid = "bounds_"+lonName+appendDomid;
884                 if (!domain->bounds_lat_name.isEmpty()) bounds_latid = domain->bounds_lat_name;
885                 else bounds_latid = "bounds_"+latName+appendDomid;
886
887                 if (domain->hasLonLat)
888                 {
889                   SuperClassWriter::addVariable(latid, typePrec, dim0, compressionLevel);
890                   SuperClassWriter::addVariable(lonid, typePrec, dim0, compressionLevel);
891                   this->writeAxisAttributes(lonid, "", "longitude", "Longitude", "degrees_east", domid);
892                   if (domain->hasBounds) SuperClassWriter::addAttribute("bounds",bounds_lonid, &lonid);
893                   this->writeAxisAttributes(latid, "", "latitude", "Latitude", "degrees_north", domid);
894                   if (domain->hasBounds) SuperClassWriter::addAttribute("bounds",bounds_latid, &latid);
895                   if (domain->hasBounds) SuperClassWriter::addDimension(dimVertId, domain->nvertex);
896                 }
897                 dim0.clear();
898                 if (domain->hasBounds)
899                 {
900                   dim0.push_back(dimXid);
901                   dim0.push_back(dimVertId);
902                   SuperClassWriter::addVariable(bounds_lonid, typePrec, dim0, compressionLevel);
903                   SuperClassWriter::addVariable(bounds_latid, typePrec, dim0, compressionLevel);
904                 }
905
906                 dim0.clear();
907                 dim0.push_back(dimXid);
908                 if (domain->hasArea)
909                 {
910                   SuperClassWriter::addVariable(areaId, typePrec, dim0, compressionLevel);
911                   SuperClassWriter::addAttribute("standard_name", StdString("cell_area"), &areaId);
912                   SuperClassWriter::addAttribute("units", StdString("m2"), &areaId);
913                 }
914
915                 SuperClassWriter::definition_end();
916
917                 if (domain->hasLonLat)
918                 {
919                   SuperClassWriter::writeData(domain->latvalue, latid, isCollective, 0);
920                   SuperClassWriter::writeData(domain->lonvalue, lonid, isCollective, 0);
921                   if (domain->hasBounds)
922                   {
923                     SuperClassWriter::writeData(domain->bounds_lonvalue, bounds_lonid, isCollective, 0);
924                     SuperClassWriter::writeData(domain->bounds_latvalue, bounds_latid, isCollective, 0);
925                   }
926                 }
927
928                 if (domain->hasArea)
929                   SuperClassWriter::writeData(domain->areavalue, areaId, isCollective, 0);
930
931                 SuperClassWriter::definition_start();
932                 break ;
933              }
934
935              case (ONE_FILE) :
936              {
937                 lonid = lonName+appendDomid;
938                 latid = latName+appendDomid;
939                 if (!domain->bounds_lon_name.isEmpty()) bounds_lonid = domain->bounds_lon_name;
940                 else bounds_lonid = "bounds_"+lonName+appendDomid;
941                 if (!domain->bounds_lat_name.isEmpty()) bounds_latid = domain->bounds_lat_name;
942                 else bounds_latid = "bounds_"+latName+appendDomid;
943
944                 dim0.push_back(dimXid);
945                 SuperClassWriter::addDimension(dimXid, domain->ni_glo);
946                 if (domain->hasLonLat)
947                 {
948                   SuperClassWriter::addVariable(latid, typePrec, dim0, compressionLevel);
949                   SuperClassWriter::addVariable(lonid, typePrec, dim0, compressionLevel);
950
951                   this->writeAxisAttributes(lonid, "", "longitude", "Longitude", "degrees_east", domid);
952                   if (domain->hasBounds) SuperClassWriter::addAttribute("bounds",bounds_lonid, &lonid);
953                   this->writeAxisAttributes(latid, "", "latitude", "Latitude", "degrees_north", domid);
954                   if (domain->hasBounds) SuperClassWriter::addAttribute("bounds",bounds_latid, &latid);
955                   if (domain->hasBounds) SuperClassWriter::addDimension(dimVertId, nvertex);
956                 }
957                 dim0.clear();
958
959                 if (domain->hasBounds)
960                 {
961                   dim0.push_back(dimXid);
962                   dim0.push_back(dimVertId);
963                   SuperClassWriter::addVariable(bounds_lonid, typePrec, dim0, compressionLevel);
964                   SuperClassWriter::addVariable(bounds_latid, typePrec, dim0, compressionLevel);
965                 }
966
967                 if (domain->hasArea)
968                 {
969                   dim0.clear();
970                   dim0.push_back(dimXid);
971                   SuperClassWriter::addVariable(areaId, typePrec, dim0, compressionLevel);
972                   SuperClassWriter::addAttribute("standard_name", StdString("cell_area"), &areaId);
973                   SuperClassWriter::addAttribute("units", StdString("m2"), &areaId);
974                 }
975
976                 SuperClassWriter::definition_end();
977
978                 std::vector<StdSize> start(1), startBounds(2) ;
979                 std::vector<StdSize> count(1), countBounds(2) ;
980                 if (domain->isEmpty())
981                 {
982                   start[0]=0 ;
983                   count[0]=0 ;
984                   startBounds[1]=0 ;
985                   countBounds[1]=nvertex ;
986                   startBounds[0]=0 ;
987                   countBounds[0]=0 ;
988                 }
989                 else
990                 {
991                   start[0]=domain->ibegin;
992                   count[0]=domain->ni;
993                   startBounds[0]=domain->ibegin;
994                   startBounds[1]=0 ;
995                   countBounds[0]=domain->ni;
996                   countBounds[1]=nvertex ;
997                 }
998
999                 if (domain->hasLonLat)
1000                 {
1001                   SuperClassWriter::writeData(domain->latvalue, latid, isCollective, 0,&start,&count);
1002                   SuperClassWriter::writeData(domain->lonvalue, lonid, isCollective, 0,&start,&count);
1003                   if (domain->hasBounds)
1004                   {
1005                     SuperClassWriter::writeData(domain->bounds_lonvalue, bounds_lonid, isCollective, 0,&startBounds,&countBounds);
1006                     SuperClassWriter::writeData(domain->bounds_latvalue, bounds_latid, isCollective, 0,&startBounds,&countBounds);
1007                   }
1008                 }
1009
1010                 if (domain->hasArea)
1011                   SuperClassWriter::writeData(domain->areavalue, areaId, isCollective, 0, &start, &count);
1012
1013                 SuperClassWriter::definition_start();
1014
1015                 break;
1016              }
1017              default :
1018                 ERROR("CNc4DataOutput::writeDomain(domain)",
1019                       << "[ type = " << SuperClass::type << "]"
1020                       << " not implemented yet !");
1021           }
1022         }
1023         catch (CNetCdfException& e)
1024         {
1025           StdString msg("On writing the domain : ");
1026           msg.append(domid); msg.append("\n");
1027           msg.append("In the context : ");
1028           msg.append(context->getId()); msg.append("\n");
1029           msg.append(e.what());
1030           ERROR("CNc4DataOutput::writeUnstructuredDomain(CDomain* domain)", << msg);
1031         }
1032         domain->addRelFile(this->filename);
1033      }
1034      //--------------------------------------------------------------
1035
1036      void CNc4DataOutput::writeAxis_(CAxis* axis)
1037      {
1038        if (axis->IsWritten(this->filename)) return;
1039
1040        // Check that the name associated to the current element is not in conflict with an existing element (due to CGrid::duplicateSentGrid)
1041        if (!axis->value.isEmpty() )
1042        {
1043          // The hash of the element will be associated to the default element name (= map key), and to the name really written
1044          int globalHash = axis->computeAttributesHash( comm_file ); // Need a MPI_Comm to distribute without redundancy some attributs (value)
1045
1046          StdString defaultNameKey = axis->getAxisOutputName();
1047          if ( !relElements_.count ( defaultNameKey ) )
1048          {
1049            // if defaultNameKey not in the map, write the element such as it is defined
1050            relElements_.insert( make_pair( defaultNameKey, make_pair(globalHash, defaultNameKey) ) );
1051          }
1052          else // look if a hash associated this key is equal
1053          {
1054            bool elementIsInMap(false);
1055            auto defaultNameKeyElements = relElements_.equal_range( defaultNameKey );
1056            for (auto it = defaultNameKeyElements.first; it != defaultNameKeyElements.second; it++)
1057            {
1058              if ( it->second.first == globalHash )
1059              {
1060                // if yes, associate the same ids to current element
1061                axis->name = it->second.second;
1062                elementIsInMap = true;
1063              }
1064            }
1065             // if no : inheritance has been excessive, define new names and store it (could be used by another grid)
1066            if (!elementIsInMap)  // ! in MAP
1067            {
1068              axis->name =  axis->getId();
1069              relElements_.insert( make_pair( defaultNameKey, make_pair(globalHash, axis->getAxisOutputName()) ) ) ;// = axis->getId()         
1070            }
1071          }
1072        }
1073
1074        axis->checkAttributes();
1075
1076        int size  = (MULTI_FILE == SuperClass::type) ? axis->n.getValue()
1077                                                          : axis->n_glo.getValue();
1078
1079        if ((0 == axis->n) && (MULTI_FILE == SuperClass::type)) return;
1080
1081        std::vector<StdString> dims;
1082        StdString axisid = axis->getAxisOutputName();
1083        StdString axisDim, axisBoundsId;
1084        if (isWrittenAxis(axisid)) return ;
1085        else setWrittenAxis(axisid);
1086
1087        nc_type typePrec ;
1088        if (axis->prec.isEmpty()) typePrec =  NC_FLOAT ;
1089        else if (axis->prec==4)   typePrec =  NC_FLOAT ;
1090        else if (axis->prec==8)   typePrec =  NC_DOUBLE ;
1091         
1092        if (!axis->label.isEmpty()) typePrec = NC_CHAR ;
1093        string strId="str_len" ;
1094        try
1095        {
1096          if (axis->dim_name.isEmpty()) axisDim = axisid;
1097          else axisDim=axis->dim_name.getValue();
1098          SuperClassWriter::addDimension(axisDim, size);
1099          dims.push_back(axisDim);
1100
1101          if (!axis->label.isEmpty() && !SuperClassWriter::dimExist(strId)) SuperClassWriter::addDimension(strId, stringArrayLen);
1102
1103          if (axis->hasValue || !axis->label.isEmpty())
1104          {
1105            if (!axis->label.isEmpty()) dims.push_back(strId);
1106
1107            SuperClassWriter::addVariable(axisid, typePrec, dims, compressionLevel);
1108
1109            if (!axis->name.isEmpty())
1110              SuperClassWriter::addAttribute("name", axis->name.getValue(), &axisid);
1111
1112            if (!axis->standard_name.isEmpty())
1113              SuperClassWriter::addAttribute("standard_name", axis->standard_name.getValue(), &axisid);
1114
1115            if (!axis->long_name.isEmpty())
1116              SuperClassWriter::addAttribute("long_name", axis->long_name.getValue(), &axisid);
1117
1118            if (!axis->unit.isEmpty())
1119              SuperClassWriter::addAttribute("units", axis->unit.getValue(), &axisid);
1120
1121            if (!axis->axis_type.isEmpty())
1122            {
1123              switch(axis->axis_type)
1124              {
1125              case CAxis::axis_type_attr::X :
1126                SuperClassWriter::addAttribute("axis", string("X"), &axisid);
1127                break;
1128              case CAxis::axis_type_attr::Y :
1129                SuperClassWriter::addAttribute("axis", string("Y"), &axisid);
1130                break;
1131              case CAxis::axis_type_attr::Z :
1132                SuperClassWriter::addAttribute("axis", string("Z"), &axisid);
1133                break;
1134              case CAxis::axis_type_attr::T :
1135                SuperClassWriter::addAttribute("axis", string("T"), &axisid);
1136                break;
1137              }
1138            }
1139
1140            if (!axis->positive.isEmpty())
1141            {
1142              SuperClassWriter::addAttribute("positive",
1143                                             (axis->positive == CAxis::positive_attr::up) ? string("up") : string("down"),
1144                                             &axisid);
1145            }
1146
1147            if (!axis->formula.isEmpty())
1148              SuperClassWriter::addAttribute("formula", axis->formula.getValue(), &axisid);
1149
1150            if (!axis->formula_term.isEmpty())
1151              SuperClassWriter::addAttribute("formula_term", axis->formula_term.getValue(), &axisid);
1152             
1153            axisBoundsId = (axis->bounds_name.isEmpty()) ? axisid + "_bounds" : axis->bounds_name;
1154            if (!axis->bounds.isEmpty() && axis->label.isEmpty())
1155            {
1156              dims.push_back("axis_nbounds");
1157              SuperClassWriter::addVariable(axisBoundsId, typePrec, dims, compressionLevel);
1158              SuperClassWriter::addAttribute("bounds", axisBoundsId, &axisid);
1159
1160              if (!axis->standard_name.isEmpty())
1161                SuperClassWriter::addAttribute("standard_name", axis->standard_name.getValue(), &axisBoundsId);
1162
1163              if (!axis->unit.isEmpty())
1164                SuperClassWriter::addAttribute("units", axis->unit.getValue(), &axisBoundsId);
1165
1166              if (!axis->formula_bounds.isEmpty())
1167                SuperClassWriter::addAttribute("formula", axis->formula_bounds.getValue(), &axisBoundsId);
1168
1169              if (!axis->formula_term_bounds.isEmpty())
1170                SuperClassWriter::addAttribute("formula_term", axis->formula_term_bounds.getValue(), &axisBoundsId);
1171            }
1172          }
1173
1174          SuperClassWriter::definition_end();
1175         
1176          switch (SuperClass::type)
1177          {
1178            case MULTI_FILE:
1179            {
1180              if (axis->label.isEmpty())
1181              {
1182                if (!axis->value.isEmpty())
1183                  SuperClassWriter::writeData(axis->value, axisid, isCollective, 0);
1184
1185                if (!axis->bounds.isEmpty())
1186                  SuperClassWriter::writeData(axis->bounds, axisBoundsId, isCollective, 0);
1187              }
1188              else
1189                SuperClassWriter::writeData(axis->label, axisid, isCollective, 0);
1190
1191              SuperClassWriter::definition_start();
1192              break;
1193            }
1194            case ONE_FILE:
1195            {
1196              std::vector<StdSize> start(1), startBounds(2) ;
1197              std::vector<StdSize> count(1), countBounds(2) ;
1198              start[0] = startBounds[0] = axis->begin;
1199              count[0] = countBounds[0] = axis->n;
1200              startBounds[1] = 0;
1201              countBounds[1] = 2;
1202
1203              if (axis->label.isEmpty())
1204              {
1205                if (!axis->value.isEmpty())
1206                  SuperClassWriter::writeData(axis->value, axisid, isCollective, 0, &start, &count);
1207
1208                if (!axis->bounds.isEmpty())
1209                  SuperClassWriter::writeData(axis->bounds, axisBoundsId, isCollective, 0, &startBounds, &countBounds);
1210              }
1211              else
1212              {
1213                std::vector<StdSize> startLabel(2), countLabel(2);
1214                startLabel[0] = start[0]; startLabel[1] = 0;
1215                countLabel[0] = count[0]; countLabel[1] = stringArrayLen;
1216                SuperClassWriter::writeData(axis->label, axisid, isCollective, 0, &startLabel, &countLabel);
1217              }
1218
1219              SuperClassWriter::definition_start();
1220
1221              break;
1222            }
1223            default :
1224              ERROR("CNc4DataOutput::writeAxis_(CAxis* axis)",
1225                    << "[ type = " << SuperClass::type << "]"
1226                    << " not implemented yet !");
1227          }
1228        }
1229        catch (CNetCdfException& e)
1230        {
1231          StdString msg("On writing the axis : ");
1232          msg.append(axisid); msg.append("\n");
1233          msg.append("In the context : ");
1234          CContext* context = CContext::getCurrent() ;
1235          msg.append(context->getId()); msg.append("\n");
1236          msg.append(e.what());
1237          ERROR("CNc4DataOutput::writeAxis_(CAxis* axis)", << msg);
1238        }
1239        axis->addRelFile(this->filename);
1240     }
1241
1242      void CNc4DataOutput::writeScalar_(CScalar* scalar)
1243      {
1244        if (scalar->IsWritten(this->filename)) return;
1245        scalar->checkAttributes();
1246        int scalarSize = 1;
1247
1248        StdString scalaId = scalar->getScalarOutputName();
1249        StdString boundsId;
1250        if (isWrittenAxis(scalaId)) return ;
1251        else setWrittenAxis(scalaId);
1252
1253        nc_type typePrec ;
1254        if (scalar->prec.isEmpty()) typePrec =  NC_FLOAT ;
1255        else if (scalar->prec==4)  typePrec =  NC_FLOAT ;
1256        else if (scalar->prec==8)   typePrec =  NC_DOUBLE ;
1257
1258        if (!scalar->label.isEmpty()) typePrec = NC_CHAR ;
1259        string strId="str_len" ;
1260
1261        try
1262        {
1263          if (!scalar->label.isEmpty() && !SuperClassWriter::dimExist(strId)) SuperClassWriter::addDimension(strId, stringArrayLen);
1264
1265          if (!scalar->value.isEmpty() || !scalar->label.isEmpty())
1266          {
1267            std::vector<StdString> dims;
1268            StdString scalarDim = scalaId;
1269
1270            if (!scalar->label.isEmpty()) dims.push_back(strId);
1271
1272            SuperClassWriter::addVariable(scalaId, typePrec, dims);
1273
1274            if (!scalar->name.isEmpty())
1275              SuperClassWriter::addAttribute("name", scalar->name.getValue(), &scalaId);
1276
1277            if (!scalar->standard_name.isEmpty())
1278              SuperClassWriter::addAttribute("standard_name", scalar->standard_name.getValue(), &scalaId);
1279
1280            if (!scalar->long_name.isEmpty())
1281              SuperClassWriter::addAttribute("long_name", scalar->long_name.getValue(), &scalaId);
1282
1283            if (!scalar->unit.isEmpty())
1284              SuperClassWriter::addAttribute("units", scalar->unit.getValue(), &scalaId);
1285
1286            if (!scalar->axis_type.isEmpty())
1287            {
1288              switch(scalar->axis_type)
1289              {
1290              case CScalar::axis_type_attr::X :
1291                SuperClassWriter::addAttribute("axis", string("X"), &scalaId);
1292                break;
1293              case CScalar::axis_type_attr::Y :
1294                SuperClassWriter::addAttribute("axis", string("Y"), &scalaId);
1295                break;
1296              case CScalar::axis_type_attr::Z :
1297                SuperClassWriter::addAttribute("axis", string("Z"), &scalaId);
1298                break;
1299              case CScalar::axis_type_attr::T :
1300                SuperClassWriter::addAttribute("axis", string("T"), &scalaId);
1301                break;
1302              }
1303            }
1304
1305            if (!scalar->positive.isEmpty())
1306            {
1307              SuperClassWriter::addAttribute("positive",
1308                                             (scalar->positive == CScalar::positive_attr::up) ? string("up") : string("down"),
1309                                             &scalaId);
1310            }
1311
1312            if (!scalar->bounds.isEmpty() && scalar->label.isEmpty())
1313            {
1314              dims.clear();
1315              dims.push_back("axis_nbounds");
1316              boundsId = (scalar->bounds_name.isEmpty()) ? (scalaId + "_bounds") : scalar->bounds_name.getValue();
1317              SuperClassWriter::addVariable(boundsId, typePrec, dims);
1318              SuperClassWriter::addAttribute("bounds", boundsId, &scalaId);
1319            }
1320
1321            SuperClassWriter::definition_end();
1322
1323            switch (SuperClass::type)
1324            {
1325              case MULTI_FILE:
1326              {
1327                CArray<double,1> scalarValue(scalarSize);
1328                CArray<string,1> scalarLabel(scalarSize);
1329                CArray<double,1> scalarBounds(scalarSize*2);
1330
1331                if (!scalar->value.isEmpty() && scalar->label.isEmpty())
1332                {
1333                  scalarValue(0) = scalar->value;
1334                  SuperClassWriter::writeData(scalarValue, scalaId, isCollective, 0);
1335                }
1336
1337                if (!scalar->bounds.isEmpty() && scalar->label.isEmpty())
1338                {
1339                  scalarBounds(0) = scalar->bounds(0);
1340                  scalarBounds(1) = scalar->bounds(1);
1341                  SuperClassWriter::writeData(scalarBounds, boundsId, isCollective, 0);
1342                }
1343
1344                if (!scalar->label.isEmpty())
1345                {
1346                  scalarLabel(0) = scalar->label;
1347                  SuperClassWriter::writeData(scalarLabel, scalaId, isCollective, 0);
1348                }
1349
1350                SuperClassWriter::definition_start();
1351
1352                break;
1353              }
1354              case ONE_FILE:
1355              {
1356                CArray<double,1> scalarValue(scalarSize);
1357                CArray<string,1> scalarLabel(scalarSize);
1358                CArray<double,1> scalarBounds(scalarSize*2);
1359
1360                std::vector<StdSize> start(1);
1361                std::vector<StdSize> count(1);
1362                start[0] = 0;
1363                count[0] = 1;
1364                if (!scalar->value.isEmpty() && scalar->label.isEmpty())
1365                {
1366                  scalarValue(0) = scalar->value;
1367                  SuperClassWriter::writeData(scalarValue, scalaId, isCollective, 0, &start, &count);
1368                }
1369                if (!scalar->bounds.isEmpty() && scalar->label.isEmpty())
1370                {
1371                  scalarBounds(0) = scalar->bounds(0);
1372                  scalarBounds(1) = scalar->bounds(1);
1373                  count[0] = 2;
1374                  SuperClassWriter::writeData(scalarBounds, boundsId, isCollective, 0, &start, &count);
1375                }
1376                if (!scalar->label.isEmpty())
1377                {
1378                  scalarLabel(0) = scalar->label;
1379                  count[0] = stringArrayLen;
1380                  SuperClassWriter::writeData(scalarLabel, scalaId, isCollective, 0, &start, &count);
1381                }
1382
1383                SuperClassWriter::definition_start();
1384
1385                break;
1386              }
1387              default :
1388                ERROR("CNc4DataOutput::writeAxis_(CAxis* scalar)",
1389                      << "[ type = " << SuperClass::type << "]"
1390                      << " not implemented yet !");
1391            }
1392          }
1393        }
1394        catch (CNetCdfException& e)
1395        {
1396          StdString msg("On writing the scalar : ");
1397          msg.append(scalaId); msg.append("\n");
1398          msg.append("In the context : ");
1399          CContext* context = CContext::getCurrent() ;
1400          msg.append(context->getId()); msg.append("\n");
1401          msg.append(e.what());
1402          ERROR("CNc4DataOutput::writeScalar_(CScalar* scalar)", << msg);
1403        }
1404        scalar->addRelFile(this->filename);
1405     }
1406
1407     //--------------------------------------------------------------
1408
1409     void CNc4DataOutput::writeGridCompressed_(CGrid* grid)
1410     {
1411        if (grid->isScalarGrid() || grid->isWrittenCompressed(this->filename)) return;
1412       
1413        // NOTA : The cuurent algorithm to write compress elements of the grid
1414        //        will work pretting well when on server side you dont't get
1415        //        partial overlap on elements between differents participating process
1416        //        So the element must be totally distributed or non distributed
1417        //        If an element is partially overlaping betwwen process then the
1418        //        total compressed part will apear artificially greater than expected
1419        //        For the current implementation of writer which is decomposed only on
1420        //        one element, it will work as expected, but for future, it must be
1421        //        reconsidered again.
1422        try
1423        {
1424          CArray<int,1> axisDomainOrder = grid->axis_domain_order;
1425          std::vector<StdString> domainList = grid->getDomainList();
1426          std::vector<StdString> axisList   = grid->getAxisList();
1427          std::vector<StdString> scalarList = grid->getScalarList();
1428          int numElement = axisDomainOrder.numElements(), idxDomain = 0, idxAxis = 0, idxScalar = 0;
1429          int commRank ;
1430          MPI_Comm_rank(comm_file,&commRank) ;
1431
1432          std::vector<StdString> dims;
1433
1434          for (int i = 0; i < numElement; ++i)
1435          {
1436            StdString varId, compress;
1437            CArray<size_t, 1> indexes;
1438            bool isDistributed;
1439            size_t nbIndexes, totalNbIndexes, offset;
1440            size_t firstGlobalIndex;
1441           
1442            if (2 == axisDomainOrder(i))
1443            {
1444              CDomain* domain = CDomain::get(domainList[idxDomain]);
1445              StdString domId = domain->getDomainOutputName();
1446
1447              if (!domain->isCompressible()
1448                  || domain->type == CDomain::type_attr::unstructured
1449                  || domain->isWrittenCompressed(this->filename)
1450                  || isWrittenCompressedDomain(domId))
1451                continue;
1452           
1453              // unstructured grid seems not be taken into account why ?
1454
1455              string lonName,latName ;
1456
1457              if (domain->lon_name.isEmpty())
1458              { 
1459                if (domain->type==CDomain::type_attr::curvilinear) lonName = "nav_lon";
1460                else lonName = "lon";
1461              }
1462              else lonName = domain->lon_name;
1463
1464              if (domain->lat_name.isEmpty())
1465              {
1466                if (domain->type==CDomain::type_attr::curvilinear) latName = "nav_lat";
1467                else latName = "lat";
1468              }
1469              else latName = domain->lat_name;
1470             
1471              StdString appendDomId  = singleDomain ? "" : "_" + domId;
1472
1473              varId = domId + "_points";
1474              compress = latName + appendDomId + " " + lonName + appendDomId;
1475     
1476              shared_ptr<CLocalView> workflowView = domain->getLocalView(CElementView::WORKFLOW) ;
1477              workflowView->getGlobalIndexView(indexes) ;
1478              nbIndexes = workflowView->getSize() ;
1479              isDistributed = domain->isDistributed();
1480              if (isDistributed)
1481              {
1482                MPI_Exscan(&nbIndexes, &offset, 1, MPI_SIZE_T, MPI_SUM, comm_file) ;
1483                if (commRank==0) offset=0 ;
1484                MPI_Allreduce(&nbIndexes,&totalNbIndexes,1 , MPI_SIZE_T, MPI_SUM, comm_file) ;
1485              }
1486              else
1487              {
1488                offset=0 ;
1489                totalNbIndexes = nbIndexes ;
1490              }
1491
1492              firstGlobalIndex = domain->ibegin + domain->jbegin * domain->ni_glo;
1493
1494              domain->addRelFileCompressed(this->filename);
1495              setWrittenCompressedDomain(domId);
1496              ++idxDomain;
1497            }
1498            else if (1 == axisDomainOrder(i))
1499            {
1500              CAxis* axis = CAxis::get(axisList[idxAxis]);
1501              StdString axisId = axis->getAxisOutputName();
1502
1503              if (!axis->isCompressible()
1504                  || axis->isWrittenCompressed(this->filename)
1505                  || isWrittenCompressedAxis(axisId))
1506                continue;
1507
1508              varId = axisId + "_points";
1509              compress = axisId;
1510
1511              shared_ptr<CLocalView> workflowView = axis->getLocalView(CElementView::WORKFLOW) ;
1512              workflowView->getGlobalIndexView(indexes) ;
1513              nbIndexes = workflowView->getSize() ;
1514              isDistributed = axis->isDistributed();
1515              if (isDistributed)
1516              {
1517                MPI_Exscan(&nbIndexes, &offset, 1, MPI_SIZE_T, MPI_SUM, comm_file) ;
1518                if (commRank==0) offset=0 ;
1519                MPI_Allreduce(&nbIndexes,&totalNbIndexes,1 , MPI_SIZE_T, MPI_SUM, comm_file) ;
1520              }
1521              else
1522              {
1523                offset=0 ;
1524                totalNbIndexes = nbIndexes ;
1525              }
1526              firstGlobalIndex = axis->begin;
1527             
1528              axis->addRelFileCompressed(this->filename);
1529              setWrittenCompressedAxis(axisId);
1530              ++idxAxis;
1531            }
1532            else
1533            {
1534              //for scalar
1535            }
1536
1537            if (!varId.empty())
1538            {
1539              SuperClassWriter::addDimension(varId, (SuperClass::type == MULTI_FILE) ? nbIndexes : totalNbIndexes);
1540
1541              dims.clear();
1542              dims.push_back(varId);
1543              SuperClassWriter::addVariable(varId, NC_UINT64, dims);
1544
1545              SuperClassWriter::addAttribute("compress", compress, &varId);
1546
1547              switch (SuperClass::type)
1548              {
1549                case (MULTI_FILE):
1550                {
1551                  indexes -= firstGlobalIndex;
1552                  SuperClassWriter::writeData(indexes, varId, isCollective, 0);
1553                  break;
1554                }
1555                case (ONE_FILE):
1556                {
1557                  std::vector<StdSize> start, count;
1558                  start.push_back(offset);
1559                  count.push_back(nbIndexes);
1560
1561                  SuperClassWriter::writeData(indexes, varId, isCollective, 0, &start, &count);
1562                  break;
1563                }
1564              }
1565            }
1566          }
1567
1568          grid->addRelFileCompressed(this->filename);
1569        }
1570        catch (CNetCdfException& e)
1571        {
1572          StdString msg("On writing compressed grid : ");
1573          msg.append(grid->getId()); msg.append("\n");
1574          msg.append("In the context : ");
1575          CContext* context = CContext::getCurrent();
1576          msg.append(context->getId()); msg.append("\n");
1577          msg.append(e.what());
1578          ERROR("CNc4DataOutput::writeGridCompressed_(CGrid* grid)", << msg);
1579        }
1580      }
1581
1582     //--------------------------------------------------------------
1583
1584     void CNc4DataOutput::writeTimeDimension_(void)
1585     {
1586       try
1587       {
1588        SuperClassWriter::addDimension(getTimeCounterName());
1589       }
1590       catch (CNetCdfException& e)
1591       {
1592         StdString msg("On writing time dimension : time_couter\n");
1593         msg.append("In the context : ");
1594         CContext* context = CContext::getCurrent() ;
1595         msg.append(context->getId()); msg.append("\n");
1596         msg.append(e.what());
1597         ERROR("CNc4DataOutput::writeTimeDimension_(void)", << msg);
1598       }
1599     }
1600
1601      //--------------------------------------------------------------
1602
1603      void CNc4DataOutput::writeField_(CField* field)
1604      {
1605        CContext* context = CContext::getCurrent() ;
1606
1607        std::vector<StdString> dims, coodinates;
1608        CGrid* grid = field->getGrid();
1609        if (!grid->doGridHaveDataToWrite())
1610          if (SuperClass::type==MULTI_FILE) return ;
1611
1612        CArray<int,1> axisDomainOrder = grid->axis_domain_order;
1613        int numElement = axisDomainOrder.numElements(), idxDomain = 0, idxAxis = 0, idxScalar = 0;
1614        std::vector<StdString> domainList = grid->getDomainList();
1615        std::vector<StdString> axisList   = grid->getAxisList();
1616        std::vector<StdString> scalarList = grid->getScalarList();       
1617
1618        StdString timeid  = getTimeCounterName();
1619        StdString dimXid,dimYid;
1620        std::deque<StdString> dimIdList, dimCoordList;
1621        bool hasArea = false;
1622        StdString cellMeasures = "area:";
1623        bool compressedOutput = !field->indexed_output.isEmpty() && field->indexed_output;
1624
1625        for (int i = 0; i < numElement; ++i)
1626        {
1627          if (2 == axisDomainOrder(i))
1628          {
1629            CDomain* domain = CDomain::get(domainList[idxDomain]);
1630            StdString domId = domain->getDomainOutputName();
1631            StdString appendDomId  = singleDomain ? "" : "_" + domId ;
1632            StdString lonName,latName ;
1633            StdString dimIname,dimJname ;
1634
1635            if (domain->lon_name.isEmpty())
1636            { 
1637              if (domain->type==CDomain::type_attr::curvilinear) lonName = "nav_lon";
1638              else lonName = "lon";
1639            }
1640            else lonName = domain->lon_name;
1641
1642            if (domain->lat_name.isEmpty())
1643            {
1644              if (domain->type==CDomain::type_attr::curvilinear) latName = "nav_lat";
1645              else latName = "lat";
1646            }
1647            else latName = domain->lat_name;
1648
1649            if (domain->dim_i_name.isEmpty())
1650            {
1651              if (domain->type==CDomain::type_attr::curvilinear) dimIname = "x";
1652              else if (domain->type==CDomain::type_attr::unstructured) dimIname = "cell";
1653              else dimIname = lonName;
1654            }
1655            else dimIname = domain->dim_i_name;
1656
1657            if (domain->dim_j_name.isEmpty())
1658            {
1659              if (domain->type==CDomain::type_attr::curvilinear) dimJname = "y";
1660              else dimJname = latName;
1661            }
1662            else dimJname = domain->dim_j_name;
1663       
1664            if (compressedOutput && domain->isCompressible() && domain->type != CDomain::type_attr::unstructured)
1665            {
1666              dimIdList.push_back(domId + "_points");
1667              field->setUseCompressedOutput();
1668            }
1669
1670            switch (domain->type)
1671            {
1672              case CDomain::type_attr::curvilinear:
1673                if (!compressedOutput || !domain->isCompressible())
1674                {
1675                  dimXid=dimIname+appendDomId;
1676                  dimYid=dimJname+appendDomId;
1677                  dimIdList.push_back(dimXid);
1678                  dimIdList.push_back(dimYid);
1679                }
1680                dimCoordList.push_back(lonName+appendDomId);
1681                dimCoordList.push_back(latName+appendDomId);
1682              break ;
1683              case CDomain::type_attr::rectilinear:
1684                if (!compressedOutput || !domain->isCompressible())
1685                {
1686                  dimXid     = dimIname+appendDomId;
1687                  dimYid     = dimJname+appendDomId;
1688                  dimIdList.push_back(dimXid);
1689                  dimIdList.push_back(dimYid);
1690                }
1691                if (lonName != dimIname)  dimCoordList.push_back(lonName+appendDomId);
1692                if (latName != dimJname)  dimCoordList.push_back(latName+appendDomId);
1693
1694              break ;
1695              case CDomain::type_attr::unstructured:
1696              {
1697                if (SuperClassWriter::useCFConvention)
1698                {
1699                  dimXid     = dimIname + appendDomId;
1700                  dimIdList.push_back(dimXid);
1701                  dimCoordList.push_back(lonName+appendDomId);
1702                  dimCoordList.push_back(latName+appendDomId);
1703                }
1704                else
1705                {
1706                  StdString domainName = domain->name;
1707                  if (domain->nvertex == 1)
1708                  {
1709                    dimXid     = "n" + domainName + "_node";
1710                    dimIdList.push_back(dimXid);
1711                    dimCoordList.push_back(StdString(domainName + "_node_x"));
1712                    dimCoordList.push_back(StdString(domainName + "_node_y"));
1713                  }
1714                  else if (domain->nvertex == 2)
1715                  {
1716                    dimXid     = "n" + domainName + "_edge";
1717                    dimIdList.push_back(dimXid);
1718                    dimCoordList.push_back(StdString(domainName + "_edge_x"));
1719                    dimCoordList.push_back(StdString(domainName + "_edge_y"));
1720                  }
1721                  else
1722                  {
1723                    dimXid     = "n" + domainName + "_face";
1724                    dimIdList.push_back(dimXid);
1725                    dimCoordList.push_back(StdString(domainName + "_face_x"));
1726                    dimCoordList.push_back(StdString(domainName + "_face_y"));
1727                  }
1728                }  // ugrid convention
1729              }  // case unstructured domain
1730            }
1731
1732            if (domain->hasArea)
1733            {
1734              hasArea = true;
1735              cellMeasures += " area" + appendDomId;
1736            }
1737            ++idxDomain;
1738          }
1739          else if (1 == axisDomainOrder(i))
1740          {
1741            CAxis* axis = CAxis::get(axisList[idxAxis]);
1742            StdString axisId = axis->getAxisOutputName();
1743            StdString axisDim;
1744
1745            if (axis->dim_name.isEmpty()) axisDim = axisId;
1746            else axisDim=axis->dim_name.getValue();
1747
1748            if (compressedOutput && axis->isCompressible())
1749            {
1750              dimIdList.push_back(axisDim + "_points");
1751              field->setUseCompressedOutput();
1752            }
1753            else
1754              dimIdList.push_back(axisDim);
1755
1756            if (axisDim != axisId) dimCoordList.push_back(axisId);
1757            ++idxAxis;
1758          }
1759          else
1760          {
1761            CScalar* scalar = CScalar::get(scalarList[idxScalar]);
1762            StdString scalarId = scalar->getScalarOutputName();
1763            if (!scalar->value.isEmpty() || !scalar->label.isEmpty())
1764              dimCoordList.push_back(scalarId);
1765            ++idxScalar;
1766          }
1767        }
1768
1769        StdString fieldid = field->getFieldOutputName();
1770
1771        nc_type type ;
1772        if (field->prec.isEmpty()) type =  NC_FLOAT ;
1773        else
1774        {
1775          if (field->prec==2) type = NC_SHORT ;
1776          else if (field->prec==4)  type =  NC_FLOAT ;
1777          else if (field->prec==8)   type =  NC_DOUBLE ;
1778        }
1779
1780        bool wtime   = !(!field->operation.isEmpty() && field->getOperationTimeType() == func::CFunctor::once);
1781
1782        if (wtime)
1783        {
1784          if (field->hasTimeInstant && hasTimeInstant) coodinates.push_back(string("time_instant"));
1785          else if (field->hasTimeCentered && hasTimeCentered)  coodinates.push_back(string("time_centered"));
1786          dims.push_back(timeid);
1787        }
1788
1789        while (!dimIdList.empty())
1790        {
1791          dims.push_back(dimIdList.back());
1792          dimIdList.pop_back();
1793        }
1794
1795        while (!dimCoordList.empty())
1796        {
1797          coodinates.push_back(dimCoordList.back());
1798          dimCoordList.pop_back();
1799        }
1800
1801        try
1802        {
1803           SuperClassWriter::addVariable(fieldid, type, dims);
1804
1805           if (!field->standard_name.isEmpty())
1806              SuperClassWriter::addAttribute
1807                 ("standard_name",  field->standard_name.getValue(), &fieldid);
1808
1809           if (!field->long_name.isEmpty())
1810              SuperClassWriter::addAttribute
1811                 ("long_name", field->long_name.getValue(), &fieldid);
1812
1813           if (!field->unit.isEmpty())
1814              SuperClassWriter::addAttribute
1815                 ("units", field->unit.getValue(), &fieldid);
1816
1817           // Ugrid field attributes "mesh" and "location"
1818           if (!SuperClassWriter::useCFConvention)
1819           {
1820            if (!domainList.empty())
1821            {
1822              CDomain* domain = CDomain::get(domainList[0]); // Suppose that we have only domain
1823              StdString mesh = domain->name;
1824              SuperClassWriter::addAttribute("mesh", mesh, &fieldid);
1825              StdString location;
1826              if (domain->nvertex == 1)
1827                location = "node";
1828              else if (domain->nvertex == 2)
1829                location = "edge";
1830              else if (domain->nvertex > 2)
1831                location = "face";
1832              SuperClassWriter::addAttribute("location", location, &fieldid);
1833            }
1834
1835           }
1836
1837           if (!field->valid_min.isEmpty())
1838              SuperClassWriter::addAttribute
1839                 ("valid_min", field->valid_min.getValue(), &fieldid);
1840
1841           if (!field->valid_max.isEmpty())
1842              SuperClassWriter::addAttribute
1843                 ("valid_max", field->valid_max.getValue(), &fieldid);
1844
1845            if (!field->scale_factor.isEmpty())
1846              SuperClassWriter::addAttribute
1847                 ("scale_factor", field->scale_factor.getValue(), &fieldid);
1848
1849             if (!field->add_offset.isEmpty())
1850              SuperClassWriter::addAttribute
1851                 ("add_offset", field->add_offset.getValue(), &fieldid);
1852
1853           SuperClassWriter::addAttribute
1854                 ("online_operation", field->operation.getValue(), &fieldid);
1855
1856          // write child variables as attributes
1857
1858
1859           bool alreadyAddCellMethod = false;
1860           StdString cellMethodsPrefix(""), cellMethodsSuffix("");
1861           if (!field->cell_methods.isEmpty())
1862           {
1863              StdString cellMethodString = field->cell_methods;
1864              if (field->cell_methods_mode.isEmpty() ||
1865                 (CField::cell_methods_mode_attr::overwrite == field->cell_methods_mode))
1866              {
1867                SuperClassWriter::addAttribute("cell_methods", cellMethodString, &fieldid);
1868                alreadyAddCellMethod = true;
1869              }
1870              else
1871              {
1872                switch (field->cell_methods_mode)
1873                {
1874                  case (CField::cell_methods_mode_attr::prefix):
1875                    cellMethodsPrefix = cellMethodString;
1876                    cellMethodsPrefix += " ";
1877                    break;
1878                  case (CField::cell_methods_mode_attr::suffix):
1879                    cellMethodsSuffix = " ";
1880                    cellMethodsSuffix += cellMethodString;
1881                    break;
1882                  case (CField::cell_methods_mode_attr::none):
1883                    break;
1884                  default:
1885                    break;
1886                }
1887              }
1888           }
1889
1890
1891           if (wtime)
1892           {
1893              CDuration freqOp = field->freq_op.getValue();
1894              freqOp.solveTimeStep(*context->calendar);
1895              StdString freqOpStr = freqOp.toStringUDUnits();
1896              SuperClassWriter::addAttribute("interval_operation", freqOpStr, &fieldid);
1897
1898              CDuration freqOut = field->getRelFile()->output_freq.getValue();
1899              freqOut.solveTimeStep(*context->calendar);
1900              SuperClassWriter::addAttribute("interval_write", freqOut.toStringUDUnits(), &fieldid);
1901
1902              StdString cellMethods(cellMethodsPrefix + "time: ");
1903              if (field->operation.getValue() == "instant") cellMethods += "point";
1904              else if (field->operation.getValue() == "average") cellMethods += "mean";
1905              else if (field->operation.getValue() == "accumulate") cellMethods += "sum";
1906              else cellMethods += field->operation;
1907              if (freqOp.resolve(*context->calendar) != freqOut.resolve(*context->calendar))
1908                cellMethods += " (interval: " + freqOpStr + ")";
1909              cellMethods += cellMethodsSuffix;
1910              if (!alreadyAddCellMethod)
1911                SuperClassWriter::addAttribute("cell_methods", cellMethods, &fieldid);
1912           }
1913
1914           if (hasArea)
1915             SuperClassWriter::addAttribute("cell_measures", cellMeasures, &fieldid);
1916
1917           if (!field->default_value.isEmpty())
1918           {
1919             double default_value = field->default_value.getValue();
1920             if (type == NC_DOUBLE)
1921             {
1922               SuperClassWriter::setDefaultValue(fieldid, &default_value);
1923             }
1924             else if (type == NC_SHORT)
1925             {
1926               short sdefault_value = (short)default_value;
1927               SuperClassWriter::setDefaultValue(fieldid, &sdefault_value);
1928             }
1929             else
1930             {
1931               float fdefault_value = (float)default_value;
1932               SuperClassWriter::setDefaultValue(fieldid, &fdefault_value);
1933             }
1934           }
1935           else
1936              SuperClassWriter::setDefaultValue(fieldid, (double*)NULL);
1937
1938            if (field->compression_level.isEmpty())
1939              field->compression_level = field->getRelFile()->compression_level.isEmpty() ? 0 : field->getRelFile()->compression_level;
1940            SuperClassWriter::setCompressionLevel(fieldid, field->compression_level);
1941
1942           {  // Ecriture des coordonnes
1943
1944              StdString coordstr; //boost::algorithm::join(coodinates, " ")
1945              std::vector<StdString>::iterator
1946                 itc = coodinates.begin(), endc = coodinates.end();
1947
1948              for (; itc!= endc; itc++)
1949              {
1950                 StdString & coord = *itc;
1951                 if (itc+1 != endc)
1952                       coordstr.append(coord).append(" ");
1953                 else  coordstr.append(coord);
1954              }
1955
1956              SuperClassWriter::addAttribute("coordinates", coordstr, &fieldid);
1957
1958           }
1959
1960           vector<CVariable*> listVars = field->getAllVariables() ;
1961           for (vector<CVariable*>::iterator it = listVars.begin() ;it != listVars.end(); it++) writeAttribute_(*it, fieldid) ;
1962
1963         }
1964         catch (CNetCdfException& e)
1965         {
1966           StdString msg("On writing field : ");
1967           msg.append(fieldid); msg.append("\n");
1968           msg.append("In the context : ");
1969           msg.append(context->getId()); msg.append("\n");
1970           msg.append(e.what());
1971           ERROR("CNc4DataOutput::writeField_(CField* field)", << msg);
1972         }
1973      } // writeField_()
1974
1975      //--------------------------------------------------------------
1976
1977      void CNc4DataOutput::writeFile_ (CFile* file)
1978      {
1979         StdString filename = file->getFileOutputName();
1980         StdString description = (!file->description.isEmpty())
1981                               ? file->description.getValue()
1982                               : StdString("Created by xios");
1983
1984         singleDomain = (file->nbDomains == 1);
1985
1986         StdString conv_str ;
1987         if (file->convention_str.isEmpty())
1988         {
1989            if (SuperClassWriter::useCFConvention) conv_str="CF-1.6" ;
1990            else conv_str="UGRID" ;
1991         }
1992         else conv_str=file->convention_str ;
1993           
1994         try
1995         {
1996           if (!appendMode) this->writeFileAttributes(filename, description,
1997                                                      conv_str,
1998                                                      StdString("An IPSL model"),
1999                                                      this->getTimeStamp());
2000
2001           if (!appendMode)
2002             SuperClassWriter::addDimension("axis_nbounds", 2);
2003         }
2004         catch (CNetCdfException& e)
2005         {
2006           StdString msg("On writing file : ");
2007           msg.append(filename); msg.append("\n");
2008           msg.append("In the context : ");
2009           CContext* context = CContext::getCurrent() ;
2010           msg.append(context->getId()); msg.append("\n");
2011           msg.append(e.what());
2012           ERROR("CNc4DataOutput::writeFile_ (CFile* file)", << msg);
2013         }
2014      }
2015
2016      void CNc4DataOutput::writeAttribute_ (CVariable* var, const string& fieldId)
2017      {
2018        StdString name = var->getVariableOutputName();
2019
2020        try
2021        {
2022          if (var->type.getValue() == CVariable::type_attr::t_int || var->type.getValue() == CVariable::type_attr::t_int32)
2023            addAttribute(name, var->getData<int>(), &fieldId);
2024          else if (var->type.getValue() == CVariable::type_attr::t_int16)
2025            addAttribute(name, var->getData<short int>(), &fieldId);
2026          else if (var->type.getValue() == CVariable::type_attr::t_float)
2027            addAttribute(name, var->getData<float>(), &fieldId);
2028          else if (var->type.getValue() == CVariable::type_attr::t_double)
2029            addAttribute(name, var->getData<double>(), &fieldId);
2030          else if (var->type.getValue() == CVariable::type_attr::t_string)
2031            addAttribute(name, var->getData<string>(), &fieldId);
2032          else
2033            ERROR("CNc4DataOutput::writeAttribute_ (CVariable* var, const string& fieldId)",
2034                  << "Unsupported variable of type " << var->type.getStringValue());
2035        }
2036       catch (CNetCdfException& e)
2037       {
2038         StdString msg("On writing attributes of variable with name : ");
2039         msg.append(name); msg.append("in the field "); msg.append(fieldId); msg.append("\n");
2040         msg.append("In the context : ");
2041         CContext* context = CContext::getCurrent() ;
2042         msg.append(context->getId()); msg.append("\n");
2043         msg.append(e.what());
2044         ERROR("CNc4DataOutput::writeAttribute_ (CVariable* var, const string& fieldId)", << msg);
2045       }
2046     }
2047
2048     void CNc4DataOutput::writeAttribute_ (CVariable* var)
2049     {
2050        StdString name = var->getVariableOutputName();
2051
2052        try
2053        {
2054          if (var->type.getValue() == CVariable::type_attr::t_int || var->type.getValue() == CVariable::type_attr::t_int32)
2055            addAttribute(name, var->getData<int>());
2056          else if (var->type.getValue() == CVariable::type_attr::t_int16)
2057            addAttribute(name, var->getData<short int>());
2058          else if (var->type.getValue() == CVariable::type_attr::t_float)
2059            addAttribute(name, var->getData<float>());
2060          else if (var->type.getValue() == CVariable::type_attr::t_double)
2061            addAttribute(name, var->getData<double>());
2062          else if (var->type.getValue() == CVariable::type_attr::t_string)
2063            addAttribute(name, var->getData<string>());
2064          else
2065            ERROR("CNc4DataOutput::writeAttribute_ (CVariable* var)",
2066                  << "Unsupported variable of type " << var->type.getStringValue());
2067        }
2068       catch (CNetCdfException& e)
2069       {
2070         StdString msg("On writing attributes of variable with name : ");
2071         msg.append(name); msg.append("\n");
2072         msg.append("In the context : ");
2073         CContext* context = CContext::getCurrent() ;
2074         msg.append(context->getId()); msg.append("\n");
2075         msg.append(e.what());
2076         ERROR("CNc4DataOutput::writeAttribute_ (CVariable* var)", << msg);
2077       }
2078     }
2079
2080      void CNc4DataOutput::syncFile_ (void)
2081      {
2082        try
2083        {
2084          SuperClassWriter::sync() ;
2085        }
2086        catch (CNetCdfException& e)
2087        {
2088         StdString msg("On synchronizing the write among processes");
2089         msg.append("In the context : ");
2090         CContext* context = CContext::getCurrent() ;
2091         msg.append(context->getId()); msg.append("\n");
2092         msg.append(e.what());
2093         ERROR("CNc4DataOutput::syncFile_ (void)", << msg);
2094        }
2095      }
2096
2097      void CNc4DataOutput::closeFile_ (void)
2098      {
2099        try
2100        {
2101          SuperClassWriter::close() ;
2102        }
2103        catch (CNetCdfException& e)
2104        {
2105         StdString msg("On closing file");
2106         msg.append("In the context : ");
2107         CContext* context = CContext::getCurrent() ;
2108         msg.append(context->getId()); msg.append("\n");
2109         msg.append(e.what());
2110         ERROR("CNc4DataOutput::syncFile_ (void)", << msg);
2111        }
2112
2113      }
2114
2115      //---------------------------------------------------------------
2116
2117      StdString CNc4DataOutput::getTimeStamp(void) const
2118      {
2119         const int buffer_size = 100;
2120         time_t rawtime;
2121         struct tm * timeinfo = NULL;
2122         char buffer [buffer_size];
2123         StdString formatStr;
2124         if (file->time_stamp_format.isEmpty()) formatStr="%Y-%b-%d %H:%M:%S %Z" ;
2125         else formatStr=file->time_stamp_format;
2126
2127//         time ( &rawtime );
2128//         timeinfo = localtime ( &rawtime );
2129         time ( &rawtime );
2130         timeinfo = gmtime ( &rawtime );
2131         strftime (buffer, buffer_size, formatStr.c_str(), timeinfo);
2132
2133         return (StdString(buffer));
2134      }
2135
2136      //---------------------------------------------------------------
2137
2138      int CNc4DataOutput::writeFieldData_ (CField*  field, const CArray<double,1>& data, const CDate& lastWrite,
2139                                           const CDate& currentWrite, int nstep)
2140      {
2141        CContext* context = CContext::getCurrent();
2142        CGrid* grid = field->getGrid();
2143       
2144        if (nstep<1) 
2145        {
2146          return nstep;
2147        }
2148       
2149        if (!grid->doGridHaveDataToWrite())
2150          if (SuperClass::type == MULTI_FILE || !isCollective)
2151          {
2152            return nstep;
2153          }
2154
2155        StdString fieldid = field->getFieldOutputName();
2156
2157        StdOStringStream oss;
2158        string timeAxisId;
2159        if (field->hasTimeInstant) timeAxisId = "time_instant";
2160        else if (field->hasTimeCentered) timeAxisId = "time_centered";
2161
2162        StdString timeBoundId = getTimeCounterName() + "_bounds";
2163
2164        StdString timeAxisBoundId;
2165        if (field->hasTimeInstant) timeAxisBoundId = "time_instant_bounds";
2166        else if (field->hasTimeCentered) timeAxisBoundId = "time_centered_bounds";
2167
2168        if (!field->wasWritten())
2169        {
2170          if (appendMode && field->getRelFile()->record_offset.isEmpty() && 
2171              field->getOperationTimeType() != func::CFunctor::once)
2172          {
2173            double factorUnit;
2174            if (!field->getRelFile()->time_units.isEmpty() && field->getRelFile()->time_units==CFile::time_units_attr::days)
2175            factorUnit=context->getCalendar()->getDayLengthInSeconds() ;
2176            else factorUnit=1 ;
2177            nstep = getRecordFromTime(currentWrite,factorUnit) + 1;
2178          }
2179
2180          field->setWritten();
2181        }
2182
2183
2184        CArray<double,1> time_data(1);
2185        CArray<double,1> time_data_bound(2);
2186        CArray<double,1> time_counter(1);
2187        CArray<double,1> time_counter_bound(2);
2188
2189        bool wtime = (field->getOperationTimeType() != func::CFunctor::once);
2190        bool wtimeCounter =false ;
2191        bool wtimeData =false ;
2192       
2193
2194        if (wtime)
2195        {
2196         
2197          if (field->hasTimeInstant)
2198          {
2199            time_data(0) = time_data_bound(1) = (Time) lastWrite;
2200            time_data_bound(0) = time_data_bound(1) = (Time) currentWrite;
2201            if (timeCounterType==instant)
2202            {
2203              time_counter(0) = time_data(0);
2204              time_counter_bound(0) = time_data_bound(0);
2205              time_counter_bound(1) = time_data_bound(1);
2206              wtimeCounter=true ;
2207            }
2208            if (hasTimeInstant) wtimeData=true ;
2209          }
2210          else if (field->hasTimeCentered)
2211          {
2212            time_data(0) = ((Time)currentWrite + (Time)lastWrite) / 2;
2213            time_data_bound(0) = (Time)lastWrite;
2214            time_data_bound(1) = (Time)currentWrite;
2215            if (timeCounterType==centered)
2216            {
2217              time_counter(0) = time_data(0) ;
2218              time_counter_bound(0) = time_data_bound(0) ;
2219              time_counter_bound(1) = time_data_bound(1) ;
2220              wtimeCounter=true ;
2221            }
2222            if (hasTimeCentered) wtimeData=true ;
2223          }
2224
2225          if (timeCounterType==record)
2226          {
2227            time_counter(0) = nstep - 1;
2228            time_counter_bound(0) = time_counter_bound(1) = nstep - 1;
2229            wtimeCounter=true ;
2230          }
2231
2232          if (!field->getRelFile()->time_units.isEmpty() && field->getRelFile()->time_units==CFile::time_units_attr::days)
2233          {
2234            double secByDay=context->getCalendar()->getDayLengthInSeconds() ;
2235            time_data/=secByDay;
2236            time_data_bound/=secByDay;
2237            time_counter/=secByDay;
2238            time_counter_bound/=secByDay;
2239          }
2240        }
2241
2242         bool isRoot = (context->intraCommRank_ == 0);
2243
2244         try
2245         {
2246           switch (SuperClass::type)
2247           {
2248              case (MULTI_FILE) :
2249              {
2250                 CTimer::get("Files : writing data").resume();
2251                 SuperClassWriter::writeData(data, fieldid, isCollective, nstep - 1);
2252                 CTimer::get("Files : writing data").suspend();
2253                 if (wtime)
2254                 {
2255                   CTimer::get("Files : writing time axis").resume();
2256                   if ( wtimeData)
2257                   {
2258                       SuperClassWriter::writeTimeAxisData(time_data, timeAxisId, isCollective, nstep - 1, isRoot);
2259                       SuperClassWriter::writeTimeAxisDataBounds(time_data_bound, timeAxisBoundId, isCollective, nstep - 1, isRoot);
2260                  }
2261                   if (wtimeCounter)
2262                   {
2263                     SuperClassWriter::writeTimeAxisData(time_counter, getTimeCounterName(), isCollective, nstep - 1,isRoot);
2264                     if (timeCounterType!=record) SuperClassWriter::writeTimeAxisDataBounds(time_counter_bound, timeBoundId, isCollective, nstep - 1, isRoot);
2265                   }
2266                   CTimer::get("Files : writing time axis").suspend();
2267                 }
2268                 break;
2269              }
2270              case (ONE_FILE) :
2271              {
2272
2273                std::vector<StdSize> start, count;
2274
2275                if (field->getUseCompressedOutput())
2276                {
2277                  CArray<int,1> axisDomainOrder = grid->axis_domain_order;
2278                  std::vector<StdString> domainList = grid->getDomainList();
2279                  std::vector<StdString> axisList   = grid->getAxisList();
2280                  int numElement = axisDomainOrder.numElements();
2281                  int idxDomain = domainList.size() - 1, idxAxis = axisList.size() - 1;
2282                  int idx = domainList.size() * 2 + axisList.size() - 1;
2283                  int commRank ;
2284
2285                  MPI_Comm_rank(comm_file,&commRank) ;
2286
2287                  start.reserve(idx+1);
2288                  count.reserve(idx+1);
2289
2290                  for (int i = numElement - 1; i >= 0; --i)
2291                  {
2292                    if (2 == axisDomainOrder(i))
2293                    {
2294                      CDomain* domain = CDomain::get(domainList[idxDomain]);
2295
2296                      if (domain->isCompressible())
2297                      {
2298                        size_t offset ;
2299                        size_t nbIndexes = domain->getLocalView(CElementView::WORKFLOW)->getSize() ;
2300                        if (domain->isDistributed())
2301                        {
2302                          MPI_Exscan(&nbIndexes, &offset, 1, MPI_SIZE_T, MPI_SUM, comm_file) ;
2303                          if (commRank==0) offset=0 ;
2304                        }
2305                        else offset=0 ;
2306
2307                        start.push_back(offset);
2308                        count.push_back(nbIndexes);
2309                        idx -= 2;
2310                      }
2311                      else
2312                      {
2313                        if ((domain->type) != CDomain::type_attr::unstructured)
2314                        {
2315                          start.push_back(domain->jbegin);
2316                          count.push_back(domain->nj);
2317                        }
2318                        --idx;
2319                        start.push_back(domain->ibegin);
2320                        count.push_back(domain->ni);
2321                        --idx;
2322                      }
2323                      --idxDomain;
2324                    }
2325                    else if (1 == axisDomainOrder(i))
2326                    {
2327                      CAxis* axis = CAxis::get(axisList[idxAxis]);
2328
2329                      if (axis->isCompressible())
2330                      {
2331                        size_t offset ;
2332                        size_t nbIndexes = axis->getLocalView(CElementView::WORKFLOW)->getSize() ;
2333                        if (axis->isDistributed())
2334                        {
2335                          MPI_Exscan(&nbIndexes, &offset, 1, MPI_SIZE_T, MPI_SUM, comm_file) ;
2336                          if (commRank==0) offset=0 ;
2337                        }
2338                        else offset=0 ;
2339
2340                        start.push_back(offset);
2341                        count.push_back(nbIndexes);
2342                      }
2343                      else
2344                      {
2345                        start.push_back(axis->begin);
2346                        count.push_back(axis->n);
2347                      }
2348                      --idxAxis;
2349                      --idx;
2350                    }
2351                  }
2352                }
2353                else
2354                {
2355                  CArray<int,1> axisDomainOrder = grid->axis_domain_order;
2356                  std::vector<StdString> domainList = grid->getDomainList();
2357                  std::vector<StdString> axisList   = grid->getAxisList();
2358                  std::vector<StdString> scalarList  = grid->getScalarList() ;
2359                  int numElement = axisDomainOrder.numElements();
2360                  int idxDomain = domainList.size() - 1, idxAxis = axisList.size() - 1;
2361                  int idx = domainList.size() * 2 + axisList.size() - 1;
2362
2363                  start.reserve(idx+1);
2364                  count.reserve(idx+1);
2365
2366                  for (int i = numElement - 1; i >= 0; --i)
2367                  {
2368                    if (2 == axisDomainOrder(i))
2369                    {
2370                      CDomain* domain = CDomain::get(domainList[idxDomain]);
2371                      if ((domain->type) != CDomain::type_attr::unstructured)
2372                      {
2373                        start.push_back(domain->jbegin);
2374                        count.push_back(domain->nj);
2375                      }
2376                      --idx ;
2377
2378                        start.push_back(domain->ibegin);
2379                        count.push_back(domain->ni);
2380                      --idx ;
2381                      --idxDomain;
2382                    }
2383                    else if (1 == axisDomainOrder(i))
2384                    {
2385                        CAxis* axis = CAxis::get(axisList[idxAxis]);
2386                        start.push_back(axis->begin);
2387                        count.push_back(axis->n);
2388                      --idx;
2389                      --idxAxis;
2390                    }
2391                    else
2392                    {
2393                      if (1 == axisDomainOrder.numElements())
2394                      {
2395                        CScalar* scalar = CScalar::get(scalarList[scalarList.size()-1]);
2396                        start.push_back(0);
2397                        count.push_back(scalar->n);
2398                      }
2399                      --idx;
2400                    }
2401                  }
2402                }
2403
2404
2405                CTimer::get("Files : writing data").resume();
2406                SuperClassWriter::writeData(data, fieldid, isCollective, nstep - 1, &start, &count);
2407                CTimer::get("Files : writing data").suspend();
2408
2409                 if (wtime)
2410                 {
2411                   CTimer::get("Files : writing time axis").resume();
2412                   if ( wtimeData)
2413                   {
2414                     SuperClassWriter::writeTimeAxisData(time_data, timeAxisId, isCollective, nstep - 1, isRoot);
2415                     SuperClassWriter::writeTimeAxisDataBounds(time_data_bound, timeAxisBoundId, isCollective, nstep - 1, isRoot);
2416                   }
2417                   if (wtimeCounter)
2418                   {
2419                     SuperClassWriter::writeTimeAxisData(time_counter, getTimeCounterName(), isCollective, nstep - 1,isRoot);
2420                     if (timeCounterType!=record) SuperClassWriter::writeTimeAxisDataBounds(time_counter_bound, timeBoundId, isCollective, nstep - 1, isRoot);
2421
2422                   }
2423                   CTimer::get("Files : writing time axis").suspend(); 
2424                 }
2425
2426                break;
2427              }
2428            }
2429            return nstep ;
2430         }
2431         catch (CNetCdfException& e)
2432         {
2433           StdString msg("On writing field data: ");
2434           msg.append(fieldid); msg.append("\n");
2435           msg.append("In the context : ");
2436           msg.append(context->getId()); msg.append("\n");
2437           msg.append(e.what());
2438           ERROR("CNc4DataOutput::writeFieldData_ (CField*  field)", << msg);
2439         }
2440      }
2441
2442      //---------------------------------------------------------------
2443
2444      void CNc4DataOutput::writeTimeAxis_
2445                  (CField*    field,
2446                   const std::shared_ptr<CCalendar> cal)
2447      {
2448         StdOStringStream oss;
2449         bool createInstantAxis=false ;
2450         bool createCenteredAxis=false ;
2451         bool createTimeCounterAxis=false ;
2452         
2453         if (field->getOperationTimeType() == func::CFunctor::once) return ;
2454
2455
2456         StdString axisId ;
2457         StdString axisBoundId;
2458         StdString timeid(getTimeCounterName());
2459         StdString timeBoundId("axis_nbounds");
2460
2461         StdString strTimeUnits ;
2462         if (!field->getRelFile()->time_units.isEmpty() && field->getRelFile()->time_units==CFile::time_units_attr::days) strTimeUnits="days since " ;
2463         else  strTimeUnits="seconds since " ;
2464 
2465         if (field->getOperationTimeType() == func::CFunctor::instant) field->hasTimeInstant = true;
2466         if (field->getOperationTimeType() == func::CFunctor::centered) field->hasTimeCentered = true;
2467
2468
2469         if (field->getRelFile()->time_counter.isEmpty())
2470         {
2471           if (timeCounterType==none) createTimeCounterAxis=true ;
2472           if (field->hasTimeCentered)
2473           {
2474             timeCounterType=centered ;
2475             if (!hasTimeCentered) createCenteredAxis=true ;
2476           }
2477           if (field->hasTimeInstant)
2478           {
2479             if (timeCounterType==none) timeCounterType=instant ;
2480             if (!hasTimeInstant) createInstantAxis=true ;
2481           }
2482         }
2483         else if (field->getRelFile()->time_counter==CFile::time_counter_attr::instant)
2484         {
2485           if (field->hasTimeCentered)
2486           {
2487             if (!hasTimeCentered) createCenteredAxis=true ;
2488           }
2489           if (field->hasTimeInstant)
2490           {
2491             if (timeCounterType==none) createTimeCounterAxis=true ;
2492             timeCounterType=instant ;
2493             if (!hasTimeInstant) createInstantAxis=true ;
2494           }
2495         }
2496         else if (field->getRelFile()->time_counter==CFile::time_counter_attr::centered)
2497         {
2498           if (field->hasTimeCentered)
2499           {
2500             if (timeCounterType==none) createTimeCounterAxis=true ;
2501             timeCounterType=centered ;
2502             if (!hasTimeCentered) createCenteredAxis=true ;
2503           }
2504           if (field->hasTimeInstant)
2505           {
2506             if (!hasTimeInstant) createInstantAxis=true ;
2507           }
2508         }
2509         else if (field->getRelFile()->time_counter==CFile::time_counter_attr::instant_exclusive)
2510         {
2511           if (field->hasTimeCentered)
2512           {
2513             if (!hasTimeCentered) createCenteredAxis=true ;
2514           }
2515           if (field->hasTimeInstant)
2516           {
2517             if (timeCounterType==none) createTimeCounterAxis=true ;
2518             timeCounterType=instant ;
2519           }
2520         }
2521         else if (field->getRelFile()->time_counter==CFile::time_counter_attr::centered_exclusive)
2522         {
2523           if (field->hasTimeCentered)
2524           {
2525             if (timeCounterType==none) createTimeCounterAxis=true ;
2526             timeCounterType=centered ;
2527           }
2528           if (field->hasTimeInstant)
2529           {
2530             if (!hasTimeInstant) createInstantAxis=true ;
2531           }
2532         }
2533         else if (field->getRelFile()->time_counter==CFile::time_counter_attr::exclusive)
2534         {
2535           if (field->hasTimeCentered)
2536           {
2537             if (timeCounterType==none) createTimeCounterAxis=true ;
2538             if (timeCounterType==instant) createInstantAxis=true ;
2539             timeCounterType=centered ;
2540           }
2541           if (field->hasTimeInstant)
2542           {
2543             if (timeCounterType==none)
2544             {
2545               createTimeCounterAxis=true ;
2546               timeCounterType=instant ;
2547             }
2548             if (timeCounterType==centered)
2549             {
2550               if (!hasTimeInstant) createInstantAxis=true ;
2551             }
2552           }
2553         }
2554         else if (field->getRelFile()->time_counter==CFile::time_counter_attr::none)
2555         {
2556           if (field->hasTimeCentered)
2557           {
2558             if (!hasTimeCentered) createCenteredAxis=true ;
2559           }
2560           if (field->hasTimeInstant)
2561           {
2562             if (!hasTimeInstant) createInstantAxis=true ;
2563           }
2564         }
2565         else if (field->getRelFile()->time_counter==CFile::time_counter_attr::record)
2566         {
2567           if (timeCounterType==none) createTimeCounterAxis=true ;
2568           timeCounterType=record ;
2569           if (field->hasTimeCentered)
2570           {
2571             if (!hasTimeCentered) createCenteredAxis=true ;
2572           }
2573           if (field->hasTimeInstant)
2574           {
2575             if (!hasTimeInstant) createInstantAxis=true ;
2576           }
2577         }
2578         
2579         if (createInstantAxis)
2580         {
2581           axisId="time_instant" ;
2582           axisBoundId="time_instant_bounds";
2583           hasTimeInstant=true ;
2584         }
2585
2586         if (createCenteredAxis)
2587         {
2588           axisId="time_centered" ;
2589           axisBoundId="time_centered_bounds";
2590           hasTimeCentered=true ;
2591         }
2592
2593         
2594         try
2595         {
2596            std::vector<StdString> dims;
2597           
2598            if (createInstantAxis || createCenteredAxis)
2599            {
2600              // Adding time_instant or time_centered
2601              dims.push_back(timeid);
2602              if (!SuperClassWriter::varExist(axisId))
2603              {
2604                SuperClassWriter::addVariable(axisId, NC_DOUBLE, dims);
2605
2606                CDate timeOrigin=cal->getTimeOrigin() ;
2607                StdOStringStream oss2;
2608                StdString strInitdate=oss2.str() ;
2609                StdString strTimeOrigin=timeOrigin.toString() ;
2610                this->writeTimeAxisAttributes(axisId, cal->getType(),strTimeUnits+strTimeOrigin,
2611                                              strTimeOrigin, axisBoundId);
2612             }
2613
2614             // Adding time_instant_bounds or time_centered_bounds variables
2615             if (!SuperClassWriter::varExist(axisBoundId))
2616             {
2617                dims.clear() ;
2618                dims.push_back(timeid);
2619                dims.push_back(timeBoundId);
2620                SuperClassWriter::addVariable(axisBoundId, NC_DOUBLE, dims);
2621             }
2622           }
2623
2624           if (createTimeCounterAxis)
2625           {
2626             // Adding time_counter
2627             axisId = getTimeCounterName();
2628             axisBoundId = getTimeCounterName() + "_bounds";
2629             dims.clear();
2630             dims.push_back(timeid);
2631             if (!SuperClassWriter::varExist(axisId))
2632             {
2633                SuperClassWriter::addVariable(axisId, NC_DOUBLE, dims);
2634                SuperClassWriter::addAttribute("axis", string("T"), &axisId);
2635
2636                if (field->getRelFile()->time_counter.isEmpty() || 
2637                   (field->getRelFile()->time_counter != CFile::time_counter_attr::record))
2638                {
2639                  CDate timeOrigin = cal->getTimeOrigin();
2640                  StdString strTimeOrigin = timeOrigin.toString();
2641
2642                  this->writeTimeAxisAttributes(axisId, cal->getType(),
2643                                                strTimeUnits+strTimeOrigin,
2644                                                strTimeOrigin, axisBoundId);
2645                }
2646             }
2647
2648             // Adding time_counter_bound dimension
2649             if (field->getRelFile()->time_counter.isEmpty() || (field->getRelFile()->time_counter != CFile::time_counter_attr::record))
2650             {
2651                if (!SuperClassWriter::varExist(axisBoundId))
2652                {
2653                  dims.clear();
2654                  dims.push_back(timeid);
2655                  dims.push_back(timeBoundId);
2656                  SuperClassWriter::addVariable(axisBoundId, NC_DOUBLE, dims);
2657                }
2658             }
2659           }
2660         }
2661         catch (CNetCdfException& e)
2662         {
2663           StdString msg("On writing time axis data: ");
2664           msg.append("In the context : ");
2665           CContext* context = CContext::getCurrent() ;
2666           msg.append(context->getId()); msg.append("\n");
2667           msg.append(e.what());
2668           ERROR("CNc4DataOutput::writeTimeAxis_ (CField*    field, \
2669                  const std::shared_ptr<CCalendar> cal)", << msg);
2670         }
2671      }
2672
2673      //---------------------------------------------------------------
2674
2675      void CNc4DataOutput::writeTimeAxisAttributes(const StdString & axis_name,
2676                                                   const StdString & calendar,
2677                                                   const StdString & units,
2678                                                   const StdString & time_origin,
2679                                                   const StdString & time_bounds,
2680                                                   const StdString & standard_name,
2681                                                   const StdString & long_name)
2682      {
2683         try
2684         {
2685           SuperClassWriter::addAttribute("standard_name", standard_name, &axis_name);
2686           SuperClassWriter::addAttribute("long_name",     long_name    , &axis_name);
2687           SuperClassWriter::addAttribute("calendar",      calendar     , &axis_name);
2688           SuperClassWriter::addAttribute("units",         units        , &axis_name);
2689           SuperClassWriter::addAttribute("time_origin",   time_origin  , &axis_name);
2690           SuperClassWriter::addAttribute("bounds",        time_bounds  , &axis_name);
2691         }
2692         catch (CNetCdfException& e)
2693         {
2694           StdString msg("On writing time axis Attribute: ");
2695           msg.append("In the context : ");
2696           CContext* context = CContext::getCurrent() ;
2697           msg.append(context->getId()); msg.append("\n");
2698           msg.append(e.what());
2699           ERROR("CNc4DataOutput::writeTimeAxisAttributes(const StdString & axis_name, \
2700                                                          const StdString & calendar,\
2701                                                          const StdString & units, \
2702                                                          const StdString & time_origin, \
2703                                                          const StdString & time_bounds, \
2704                                                          const StdString & standard_name, \
2705                                                          const StdString & long_name)", << msg);
2706         }
2707      }
2708
2709      //---------------------------------------------------------------
2710
2711      void CNc4DataOutput::writeAxisAttributes(const StdString & axis_name,
2712                                               const StdString & axis,
2713                                               const StdString & standard_name,
2714                                               const StdString & long_name,
2715                                               const StdString & units,
2716                                               const StdString & nav_model)
2717      {
2718         try
2719         {
2720          if (!axis.empty())
2721            SuperClassWriter::addAttribute("axis"       , axis         , &axis_name);
2722
2723          SuperClassWriter::addAttribute("standard_name", standard_name, &axis_name);
2724          SuperClassWriter::addAttribute("long_name"    , long_name    , &axis_name);
2725          SuperClassWriter::addAttribute("units"        , units        , &axis_name);
2726//          SuperClassWriter::addAttribute("nav_model"    , nav_model    , &axis_name);
2727         }
2728         catch (CNetCdfException& e)
2729         {
2730           StdString msg("On writing Axis Attribute: ");
2731           msg.append("In the context : ");
2732           CContext* context = CContext::getCurrent() ;
2733           msg.append(context->getId()); msg.append("\n");
2734           msg.append(e.what());
2735           ERROR("CNc4DataOutput::writeAxisAttributes(const StdString & axis_name, \
2736                                                      const StdString & axis, \
2737                                                      const StdString & standard_name, \
2738                                                      const StdString & long_name, \
2739                                                      const StdString & units, \
2740                                                      const StdString & nav_model)", << msg);
2741         }
2742      }
2743
2744      //---------------------------------------------------------------
2745
2746      void CNc4DataOutput::writeLocalAttributes
2747         (int ibegin, int ni, int jbegin, int nj, StdString domid)
2748      {
2749        try
2750        {
2751         SuperClassWriter::addAttribute(StdString("ibegin").append(domid), ibegin);
2752         SuperClassWriter::addAttribute(StdString("ni"    ).append(domid), ni);
2753         SuperClassWriter::addAttribute(StdString("jbegin").append(domid), jbegin);
2754         SuperClassWriter::addAttribute(StdString("nj"    ).append(domid), nj);
2755        }
2756        catch (CNetCdfException& e)
2757        {
2758           StdString msg("On writing Local Attributes: ");
2759           msg.append("In the context : ");
2760           CContext* context = CContext::getCurrent() ;
2761           msg.append(context->getId()); msg.append("\n");
2762           msg.append(e.what());
2763           ERROR("CNc4DataOutput::writeLocalAttributes \
2764                  (int ibegin, int ni, int jbegin, int nj, StdString domid)", << msg);
2765        }
2766
2767      }
2768
2769      void CNc4DataOutput::writeLocalAttributes_IOIPSL(const StdString& dimXid, const StdString& dimYid,
2770                                                       int ibegin, int ni, int jbegin, int nj, int ni_glo, int nj_glo, int rank, int size)
2771      {
2772         CArray<int,1> array(2) ;
2773
2774         try
2775         {
2776           SuperClassWriter::addAttribute("DOMAIN_number_total",size ) ;
2777           SuperClassWriter::addAttribute("DOMAIN_number", rank) ;
2778           array = SuperClassWriter::getDimension(dimXid) + 1, SuperClassWriter::getDimension(dimYid) + 1;
2779           SuperClassWriter::addAttribute("DOMAIN_dimensions_ids",array) ;
2780           array=ni_glo,nj_glo ;
2781           SuperClassWriter::addAttribute("DOMAIN_size_global", array) ;
2782           array=ni,nj ;
2783           SuperClassWriter::addAttribute("DOMAIN_size_local", array) ;
2784           array=ibegin+1,jbegin+1 ;
2785           SuperClassWriter::addAttribute("DOMAIN_position_first", array) ;
2786           array=ibegin+ni-1+1,jbegin+nj-1+1 ;
2787           SuperClassWriter::addAttribute("DOMAIN_position_last",array) ;
2788           array=0,0 ;
2789           SuperClassWriter::addAttribute("DOMAIN_halo_size_start", array) ;
2790           SuperClassWriter::addAttribute("DOMAIN_halo_size_end", array);
2791           SuperClassWriter::addAttribute("DOMAIN_type",string("box")) ;
2792  /*
2793           SuperClassWriter::addAttribute("DOMAIN_DIM_N001",string("x")) ;
2794           SuperClassWriter::addAttribute("DOMAIN_DIM_N002",string("y")) ;
2795           SuperClassWriter::addAttribute("DOMAIN_DIM_N003",string("axis_A")) ;
2796           SuperClassWriter::addAttribute("DOMAIN_DIM_N004",string("time_counter")) ;
2797  */
2798         }
2799         catch (CNetCdfException& e)
2800         {
2801           StdString msg("On writing Local Attributes IOIPSL \n");
2802           msg.append("In the context : ");
2803           CContext* context = CContext::getCurrent() ;
2804           msg.append(context->getId()); msg.append("\n");
2805           msg.append(e.what());
2806           ERROR("CNc4DataOutput::writeLocalAttributes_IOIPSL \
2807                  (int ibegin, int ni, int jbegin, int nj, int ni_glo, int nj_glo, int rank, int size)", << msg);
2808         }
2809      }
2810      //---------------------------------------------------------------
2811
2812      void CNc4DataOutput:: writeFileAttributes(const StdString & name,
2813                                                const StdString & description,
2814                                                const StdString & conventions,
2815                                                const StdString & production,
2816                                                const StdString & timeStamp)
2817      {
2818         try
2819         {
2820           SuperClassWriter::addAttribute("name"       , name);
2821           SuperClassWriter::addAttribute("description", description);
2822           SuperClassWriter::addAttribute("title"      , description);
2823           SuperClassWriter::addAttribute("Conventions", conventions);
2824           // SuperClassWriter::addAttribute("production" , production);
2825
2826           StdString timeStampStr ;
2827           if (file->time_stamp_name.isEmpty()) timeStampStr="timeStamp" ;
2828           else timeStampStr=file->time_stamp_name ;
2829           SuperClassWriter::addAttribute(timeStampStr, timeStamp);
2830
2831           StdString uuidName ;
2832           if (file->uuid_name.isEmpty()) uuidName="uuid" ;
2833           else uuidName=file->uuid_name ;
2834
2835           if (file->uuid_format.isEmpty()) SuperClassWriter::addAttribute(uuidName, getUuidStr());
2836           else SuperClassWriter::addAttribute(uuidName, getUuidStr(file->uuid_format));
2837         
2838         }
2839         catch (CNetCdfException& e)
2840         {
2841           StdString msg("On writing File Attributes \n ");
2842           msg.append("In the context : ");
2843           CContext* context = CContext::getCurrent() ;
2844           msg.append(context->getId()); msg.append("\n");
2845           msg.append(e.what());
2846           ERROR("CNc4DataOutput:: writeFileAttributes(const StdString & name, \
2847                                                const StdString & description, \
2848                                                const StdString & conventions, \
2849                                                const StdString & production, \
2850                                                const StdString & timeStamp)", << msg);
2851         }
2852      }
2853
2854      //---------------------------------------------------------------
2855
2856      void CNc4DataOutput::writeMaskAttributes(const StdString & mask_name,
2857                                               int data_dim,
2858                                               int data_ni,
2859                                               int data_nj,
2860                                               int data_ibegin,
2861                                               int data_jbegin)
2862      {
2863         try
2864         {
2865           SuperClassWriter::addAttribute("data_dim"   , data_dim   , &mask_name);
2866           SuperClassWriter::addAttribute("data_ni"    , data_ni    , &mask_name);
2867           SuperClassWriter::addAttribute("data_nj"    , data_nj    , &mask_name);
2868           SuperClassWriter::addAttribute("data_ibegin", data_ibegin, &mask_name);
2869           SuperClassWriter::addAttribute("data_jbegin", data_jbegin, &mask_name);
2870         }
2871         catch (CNetCdfException& e)
2872         {
2873           StdString msg("On writing Mask Attributes \n ");
2874           msg.append("In the context : ");
2875           CContext* context = CContext::getCurrent() ;
2876           msg.append(context->getId()); msg.append("\n");
2877           msg.append(e.what());
2878           ERROR("CNc4DataOutput::writeMaskAttributes(const StdString & mask_name, \
2879                                               int data_dim, \
2880                                               int data_ni, \
2881                                               int data_nj, \
2882                                               int data_ibegin, \
2883                                               int data_jbegin)", << msg);
2884         }
2885      }
2886
2887      ///--------------------------------------------------------------
2888
2889      StdSize CNc4DataOutput::getRecordFromTime(Time time, double factorUnit)
2890      {
2891        std::map<Time, StdSize>::const_iterator it = timeToRecordCache.find(time);
2892        if (it == timeToRecordCache.end())
2893        {
2894          StdString timeAxisBoundsId(getTimeCounterName() + "_bounds");
2895          if (!SuperClassWriter::varExist(timeAxisBoundsId)) timeAxisBoundsId = "time_centered_bounds";
2896          if (!SuperClassWriter::varExist(timeAxisBoundsId)) timeAxisBoundsId = "time_instant_bounds";
2897
2898          CArray<double,2> timeAxisBounds;
2899          std::vector<StdSize> dimSize(SuperClassWriter::getDimensions(timeAxisBoundsId)) ;
2900         
2901          StdSize record = 0;
2902          double dtime(time);
2903          for (int n = dimSize[0] - 1; n >= 0; n--)
2904          {
2905            SuperClassWriter::getTimeAxisBounds(timeAxisBounds, timeAxisBoundsId, isCollective, n);
2906            timeAxisBounds*=factorUnit ;
2907            if (timeAxisBounds(1, 0) < dtime)
2908            {
2909              record = n + 1;
2910              break;
2911            }
2912          }
2913          it = timeToRecordCache.insert(std::make_pair(time, record)).first;
2914        }
2915        return it->second;
2916      }
2917
2918      ///--------------------------------------------------------------
2919
2920      bool CNc4DataOutput::isWrittenDomain(const std::string& domainName) const
2921      {
2922        return (this->writtenDomains.find(domainName) != this->writtenDomains.end());
2923      }
2924
2925      bool CNc4DataOutput::isWrittenCompressedDomain(const std::string& domainName) const
2926      {
2927        return (this->writtenCompressedDomains.find(domainName) != this->writtenCompressedDomains.end());
2928      }
2929
2930      bool CNc4DataOutput::isWrittenAxis(const std::string& axisName) const
2931      {
2932        return (this->writtenAxis.find(axisName) != this->writtenAxis.end());
2933      }
2934
2935      bool CNc4DataOutput::isWrittenCompressedAxis(const std::string& axisName) const
2936      {
2937        return (this->writtenCompressedAxis.find(axisName) != this->writtenCompressedAxis.end());
2938      }
2939
2940      bool CNc4DataOutput::isWrittenScalar(const std::string& scalarName) const
2941      {
2942        return (this->writtenScalar.find(scalarName) != this->writtenScalar.end());
2943      }
2944
2945      void CNc4DataOutput::setWrittenDomain(const std::string& domainName)
2946      {
2947        this->writtenDomains.insert(domainName);
2948      }
2949
2950      void CNc4DataOutput::setWrittenCompressedDomain(const std::string& domainName)
2951      {
2952        this->writtenCompressedDomains.insert(domainName);
2953      }
2954
2955      void CNc4DataOutput::setWrittenAxis(const std::string& axisName)
2956      {
2957        this->writtenAxis.insert(axisName);
2958      }
2959
2960      void CNc4DataOutput::setWrittenCompressedAxis(const std::string& axisName)
2961      {
2962        this->writtenCompressedAxis.insert(axisName);
2963      }
2964
2965      void CNc4DataOutput::setWrittenScalar(const std::string& scalarName)
2966      {
2967        this->writtenScalar.insert(scalarName);
2968      }
2969} // namespace xios
Note: See TracBrowser for help on using the repository browser.