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

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

Extend domain name management using hash table to unstructured domain.

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