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

Last change on this file was 2629, checked in by jderouillat, 2 months ago

Delete boost dependencies, the few features used are replaced by functions stored in extern/boost_extraction

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