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

Last change on this file since 2391 was 2391, checked in by ymipsl, 22 months ago

Solvue issue when trying to output scalar in nectdf file coming from a grid composed of multiple scalar tensor product.
YM

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