Salome HOME
Merge remote branch 'origin/V7_dev'
[modules/paravis.git] / src / Plugins / MEDReader / IO / vtkMEDReader.cxx
1 // Copyright (C) 2010-2016  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19 // Author : Anthony Geay
20
21 #include "vtkMEDReader.h"
22 #include "vtkGenerateVectors.h"
23
24 #include "vtkMultiBlockDataSet.h"
25 #include "vtkInformation.h"
26 #include "vtkDataSetAttributes.h"
27 #include "vtkStringArray.h"
28 #include "vtkMutableDirectedGraph.h"
29 #include "vtkInformationStringKey.h"
30 //
31 #include "vtkUnsignedCharArray.h"
32 #include "vtkInformationVector.h"
33 #include "vtkSmartPointer.h"
34 #include "vtkVariantArray.h"
35 #include "vtkExecutive.h"
36 #include "vtkStreamingDemandDrivenPipeline.h"
37 #include "vtkMultiTimeStepAlgorithm.h"
38 #include "vtkUnstructuredGrid.h"
39 #include "vtkInformationQuadratureSchemeDefinitionVectorKey.h"
40 #include "vtkQuadratureSchemeDefinition.h"
41 #include "vtkPointData.h"
42 #include "vtkCellData.h"
43 #include "vtkCellType.h"
44 #include "vtkCellArray.h"
45 #include "vtkDoubleArray.h"
46 #include "vtkObjectFactory.h"
47 #include "vtkInformationDataObjectMetaDataKey.h"
48
49 #ifdef MEDREADER_USE_MPI
50 #include "vtkMultiProcessController.h"
51 #endif
52
53 #include "MEDFileFieldRepresentationTree.hxx"
54
55 #include <map>
56 #include <string>
57 #include <vector>
58 #include <sstream>
59 #include <algorithm>
60
61 /*!
62  * This class stores properties in loading state mode (pvsm) when the MED file has not been read yet.
63  * The file is not read beacause FileName has not been informed yet ! So this class stores properties of vtkMEDReader instance that
64  * owns it and wait the vtkMEDReader::SetFileName to apply properties afterwards.
65  */
66 class PropertyKeeper
67 {
68 public:
69   PropertyKeeper(vtkMEDReader *master):_master(master),IsGVActivated(false),GVValue(0),IsCMActivated(false),CMValue(0) { }
70   void assignPropertiesIfNeeded();
71   bool arePropertiesOnTreeToSetAfter() const;
72   //
73   void pushFieldStatusEntry(const char* name, int status);
74   void pushGenerateVectorsValue(int value);
75   void pushChangeModeValue(int value);
76   void pushTimesFlagsStatusEntry(const char* name, int status);
77 protected:
78   // pool of pairs to assign in SetFieldsStatus if needed. The use case is the load using pvsm.
79   std::vector< std::pair<std::string,int> > SetFieldsStatusPairs;
80   // generate vector
81   bool IsGVActivated;
82   int GVValue;
83   // change mode
84   bool IsCMActivated;
85   int CMValue;
86   //
87   std::vector< std::pair<std::string,int> > TimesFlagsStatusPairs;
88   vtkMEDReader *_master;
89 };
90
91 void PropertyKeeper::assignPropertiesIfNeeded()
92 {
93   if(!this->SetFieldsStatusPairs.empty())
94     {
95       for(std::vector< std::pair<std::string,int> >::const_iterator it=this->SetFieldsStatusPairs.begin();it!=this->SetFieldsStatusPairs.end();it++)
96         _master->SetFieldsStatus((*it).first.c_str(),(*it).second);
97       this->SetFieldsStatusPairs.clear();
98     }
99   if(!this->TimesFlagsStatusPairs.empty())
100     {
101       for(std::vector< std::pair<std::string,int> >::const_iterator it=this->TimesFlagsStatusPairs.begin();it!=this->TimesFlagsStatusPairs.end();it++)
102         _master->SetTimesFlagsStatus((*it).first.c_str(),(*it).second);
103       this->TimesFlagsStatusPairs.clear();
104     }
105   if(this->IsGVActivated)
106     {
107       _master->GenerateVectors(this->GVValue);
108       this->IsGVActivated=false;
109     }
110   if(this->IsCMActivated)
111     {
112       _master->ChangeMode(this->CMValue);
113       this->IsCMActivated=false;
114     }
115 }
116
117 void PropertyKeeper::pushFieldStatusEntry(const char* name, int status)
118 {
119   bool found(false);
120   for(std::vector< std::pair<std::string,int> >::const_iterator it=this->SetFieldsStatusPairs.begin();it!=this->SetFieldsStatusPairs.end() && !found;it++)
121     found=(*it).first==name;
122   if(!found)
123     this->SetFieldsStatusPairs.push_back(std::pair<std::string,int>(name,status));
124 }
125
126 void PropertyKeeper::pushTimesFlagsStatusEntry(const char* name, int status)
127 {
128   bool found(false);
129   for(std::vector< std::pair<std::string,int> >::const_iterator it=this->TimesFlagsStatusPairs.begin();it!=this->TimesFlagsStatusPairs.end() && !found;it++)
130     found=(*it).first==name;
131   if(!found)
132     this->TimesFlagsStatusPairs.push_back(std::pair<std::string,int>(name,status));
133 }
134
135 void PropertyKeeper::pushGenerateVectorsValue(int value)
136 {
137   this->IsGVActivated=true;
138   this->GVValue=value;
139 }
140
141 void PropertyKeeper::pushChangeModeValue(int value)
142 {
143   this->IsCMActivated=true;
144   this->CMValue=value;
145 }
146
147 bool PropertyKeeper::arePropertiesOnTreeToSetAfter() const
148 {
149   return !SetFieldsStatusPairs.empty();
150 }
151
152 class vtkMEDReader::vtkMEDReaderInternal
153 {
154
155 public:
156   vtkMEDReaderInternal(vtkMEDReader *master):TK(0),IsMEDOrSauv(true),IsStdOrMode(false),GenerateVect(false),SIL(0),LastLev0(-1),FirstCall0(2),PK(master),MyMTime(0)
157   {
158   }
159
160   bool PluginStart0()
161   {
162     return false; // TODO Useless and buggy
163     if(FirstCall0==0)
164       return false;
165     FirstCall0--;
166     return true;
167   }
168
169   ~vtkMEDReaderInternal()
170   {
171     if(this->SIL)
172       this->SIL->Delete();
173   }
174 public:
175   MEDFileFieldRepresentationTree Tree;
176   TimeKeeper TK;
177   std::string FileName;
178   //when true the file is MED file. when false it is a Sauv file
179   bool IsMEDOrSauv;
180   //when false -> std, true -> mode. By default std (false).
181   bool IsStdOrMode;
182   //when false -> do nothing. When true cut off or extend to nbOfCompo=3 vector arrays.
183   bool GenerateVect;
184   std::string DftMeshName;
185   // Store the vtkMutableDirectedGraph that represents links between family, groups and cell types
186   vtkMutableDirectedGraph* SIL;
187   // store the lev0 id in Tree corresponding to the TIME_STEPS in the pipeline.
188   int LastLev0;
189   // The property keeper is usable only in pvsm mode.
190   PropertyKeeper PK;
191   int MyMTime;
192   std::set<std::string> _wonderful_set;// this set is used by SetFieldsStatus method to detect the fact that SetFieldsStatus has been called for all items ! Great Items are not sorted ! Why ?
193   std::map<std::string,bool> _wonderful_ref;// this map stores the state before a SetFieldsStatus status.
194
195 private:
196   unsigned char FirstCall0;
197 };
198
199 vtkStandardNewMacro(vtkMEDReader);
200
201 // vtkInformationKeyMacro(vtkMEDReader, META_DATA, DataObjectMetaData); // Here we need to customize vtkMEDReader::META_DATA method
202 // start of overload of vtkInformationKeyMacro
203 static vtkInformationDataObjectMetaDataKey *vtkMEDReader_META_DATA=new vtkInformationDataObjectMetaDataKey("META_DATA","vtkMEDReader");
204
205 vtkInformationDataObjectMetaDataKey *vtkMEDReader::META_DATA()  
206 {
207   static const char ZE_KEY[]="vtkMEDReader::META_DATA";
208   vtkInformationDataObjectMetaDataKey *ret(vtkMEDReader_META_DATA);
209   MEDCoupling::GlobalDict *gd(MEDCoupling::GlobalDict::GetInstance());
210   if(!gd->hasKey(ZE_KEY))
211     {// here META_DATA is put on global var to be exchanged with other filters without dependancy of MEDReader. Please do not change ZE_KEY !
212       std::ostringstream oss; oss << ret;
213       gd->setKeyValue(ZE_KEY,oss.str());
214     }
215   return ret;
216 }
217 // end of overload of vtkInformationKeyMacro
218
219 vtkMEDReader::vtkMEDReader():Internal(new vtkMEDReaderInternal(this))
220 {
221   this->SetNumberOfInputPorts(0);
222   this->SetNumberOfOutputPorts(1);
223 }
224
225 vtkMEDReader::~vtkMEDReader()
226 {
227   delete this->Internal;
228   this->Internal = 0;
229 }
230
231 void vtkMEDReader::Reload()
232 {
233   std::string fName((const char *)this->GetFileName());
234   delete this->Internal;
235   this->Internal=new vtkMEDReaderInternal(this);
236   this->SetFileName(fName.c_str());
237 }
238
239 int vtkMEDReader::GetServerModifTime()
240 {
241   if( !this->Internal )
242     return -1;
243   return this->Internal->MyMTime;
244 }
245
246 void vtkMEDReader::GenerateVectors(int val)
247 {
248   if ( !this->Internal )
249     return;
250   
251   if(this->Internal->FileName.empty())
252     {//pvsm mode
253       this->Internal->PK.pushGenerateVectorsValue(val);
254       return ;
255     }
256   //not pvsm mode (general case)
257   bool val2((bool)val);
258   if(val2!=this->Internal->GenerateVect)
259     {
260       this->Internal->GenerateVect=val2;
261       this->Modified();
262     }
263 }
264
265 void vtkMEDReader::ChangeMode(int newMode)
266 {
267   if ( !this->Internal )
268     return;
269   
270   if(this->Internal->FileName.empty())
271     {//pvsm mode
272       this->Internal->PK.pushChangeModeValue(newMode);
273       return ;
274     }
275   //not pvsm mode (general case)
276   this->Internal->IsStdOrMode=newMode!=0;
277   this->Modified();
278 }
279
280 const char *vtkMEDReader::GetSeparator()
281 {
282   return MEDFileFieldRepresentationLeavesArrays::ZE_SEP;
283 }
284
285 void vtkMEDReader::SetFileName(const char *fname)
286 {
287   if(!this->Internal)
288     return;
289   try
290     {
291       this->Internal->FileName=fname;
292       std::size_t pos(this->Internal->FileName.find_last_of('.'));
293       if(pos!=std::string::npos)
294         {
295           std::string ext(this->Internal->FileName.substr(pos));
296           if(ext.find("sauv")!=std::string::npos)
297             this->Internal->IsMEDOrSauv=false;
298         }
299       if(this->Internal->Tree.getNumberOfLeavesArrays()==0)
300         {
301           int iPart(-1),nbOfParts(-1);
302 #ifdef MEDREADER_USE_MPI
303           vtkMultiProcessController *vmpc(vtkMultiProcessController::GetGlobalController());
304           if(vmpc)
305             {
306               iPart=vmpc->GetLocalProcessId();
307               nbOfParts=vmpc->GetNumberOfProcesses();
308             }
309 #endif
310           this->Internal->Tree.loadMainStructureOfFile(this->Internal->FileName.c_str(),this->Internal->IsMEDOrSauv,iPart,nbOfParts);
311           if(!this->Internal->PK.arePropertiesOnTreeToSetAfter())
312             this->Internal->Tree.activateTheFirst();//This line manually initialize the status of server (this) with the remote client.
313           this->Internal->TK.setMaxNumberOfTimeSteps(this->Internal->Tree.getMaxNumberOfTimeSteps());
314         }
315       this->Modified();
316       this->Internal->PK.assignPropertiesIfNeeded();
317     }
318   catch(INTERP_KERNEL::Exception& e)
319     {
320       delete this->Internal;
321       this->Internal=0;
322       std::ostringstream oss;
323       oss << "Exception has been thrown in vtkMEDReader::SetFileName : " << e.what() << std::endl;
324       if(this->HasObserver("ErrorEvent") )
325         this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
326       else
327         vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
328       vtkObject::BreakOnError();
329     }
330 }
331
332 char *vtkMEDReader::GetFileName()
333 {
334   if (!this->Internal)
335     return 0;
336   return const_cast<char *>(this->Internal->FileName.c_str());
337 }
338
339 int vtkMEDReader::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
340 {
341 //  std::cout << "########################################## vtkMEDReader::RequestInformation ##########################################" << std::endl;
342   if(!this->Internal)
343     return 0;
344   try
345     {
346 //      request->Print(cout);
347       vtkInformation *outInfo(outputVector->GetInformationObject(0));
348       outInfo->Set(vtkDataObject::DATA_TYPE_NAME(),"vtkMultiBlockDataSet");
349       this->UpdateSIL(request, outInfo);
350
351       // Set the meta data graph as a meta data key in the information
352       // That's all that is needed to transfer it along the pipeline
353       outInfo->Set(vtkMEDReader::META_DATA(),this->Internal->SIL);
354
355       bool dummy(false);
356       this->PublishTimeStepsIfNeeded(outInfo,dummy);
357     }
358   catch(INTERP_KERNEL::Exception& e)
359     {
360       std::ostringstream oss;
361       oss << "Exception has been thrown in vtkMEDReader::RequestInformation : " << e.what() << std::endl;
362       if(this->HasObserver("ErrorEvent") )
363         this->InvokeEvent("ErrorEvent",const_cast<char *>(oss.str().c_str()));
364       else
365         vtkOutputWindowDisplayErrorText(const_cast<char *>(oss.str().c_str()));
366       vtkObject::BreakOnError();
367       return 0;
368     }
369   return 1;
370 }
371
372 int vtkMEDReader::RequestData(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
373 {
374 //  std::cout << "########################################## vtkMEDReader::RequestData ##########################################" << std::endl;
375   if(!this->Internal)
376     return 0;
377   try
378     {
379 //      request->Print(cout);
380       vtkInformation *outInfo(outputVector->GetInformationObject(0));
381       vtkMultiBlockDataSet *output(vtkMultiBlockDataSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
382       bool isUpdated(false);
383       double reqTS(0.);
384       if(outInfo->Has(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP()))
385         reqTS=outInfo->Get(vtkStreamingDemandDrivenPipeline::UPDATE_TIME_STEP());
386       this->FillMultiBlockDataSetInstance(output,reqTS);
387       output->GetInformation()->Set(vtkDataObject::DATA_TIME_STEP(),reqTS);
388
389       // Is it really needed ? TODO
390       this->UpdateSIL(request, outInfo);
391     }
392   catch(INTERP_KERNEL::Exception& e)
393     {
394       std::cerr << "Exception has been thrown in vtkMEDReader::RequestData : " << e.what() << std::endl;
395       return 0;
396     }
397   return 1;
398 }
399
400 void vtkMEDReader::SetFieldsStatus(const char* name, int status)
401 {
402   if( !this->Internal )
403     return;
404
405   //this->Internal->_wonderful_set.insert(name);
406   if(this->Internal->FileName.empty())
407     {//pvsm mode
408       this->Internal->PK.pushFieldStatusEntry(name,status);
409       return ;
410     }
411   if(this->Internal->_wonderful_set.empty())
412     this->Internal->_wonderful_ref=this->Internal->Tree.dumpState();// start of SetFieldsStatus serie -> store ref to compare at the end of the SetFieldsStatus serie.
413   this->Internal->_wonderful_set.insert(name);
414   //not pvsm mode (general case)
415   try
416     {
417       this->Internal->Tree.changeStatusOfAndUpdateToHaveCoherentVTKDataSet(this->Internal->Tree.getIdHavingZeName(name),status);
418       if(this->Internal->_wonderful_set.size()==GetNumberOfFieldsTreeArrays())
419         {
420           if(this->Internal->_wonderful_ref!=this->Internal->Tree.dumpState())
421             {
422               if(!this->Internal->PluginStart0())
423                 {
424                   this->Modified();
425                 }
426               this->Internal->MyMTime++;
427             }
428           this->Internal->_wonderful_set.clear();
429         }
430     }
431   catch(INTERP_KERNEL::Exception& e)
432     {
433       if(!this->Internal->FileName.empty())
434         {
435           std::cerr << "vtkMEDReader::SetFieldsStatus error : " << e.what() << " *** WITH STATUS=" << status << endl;
436           return ;
437         }
438     }
439 }
440
441 int vtkMEDReader::GetNumberOfFieldsTreeArrays()
442 {
443   if(!this->Internal)
444     return 0;
445   return this->Internal->Tree.getNumberOfLeavesArrays();
446 }
447
448 const char *vtkMEDReader::GetFieldsTreeArrayName(int index)
449 {
450   if(!this->Internal)
451     return 0;
452   return this->Internal->Tree.getNameOfC(index);
453 }
454
455 int vtkMEDReader::GetFieldsTreeArrayStatus(const char *name)
456 {
457   if(!this->Internal)
458     return -1;
459
460   int zeId(this->Internal->Tree.getIdHavingZeName(name));
461   int ret(this->Internal->Tree.getStatusOf(zeId));
462   return ret;
463 }
464
465 void vtkMEDReader::SetTimesFlagsStatus(const char *name, int status)
466 {
467   if (!this->Internal)
468     return;
469   
470   if(this->Internal->FileName.empty())
471     {//pvsm mode
472       this->Internal->PK.pushTimesFlagsStatusEntry(name,status);
473       return ;
474     }
475   //not pvsm mode (general case)
476   int pos(0);
477   std::istringstream iss(name); iss >> pos;
478   this->Internal->TK.getTimesFlagArray()[pos].first=(bool)status;
479   if(pos==this->Internal->TK.getTimesFlagArray().size()-1)
480     if(!this->Internal->PluginStart0())
481       {
482         this->Modified();
483         //this->Internal->TK.printSelf(std::cerr);
484       }
485 }
486
487 int vtkMEDReader::GetNumberOfTimesFlagsArrays()
488 {
489   if(!this->Internal)
490     return 0;
491   return (int)this->Internal->TK.getTimesFlagArray().size();
492 }
493
494 const char *vtkMEDReader::GetTimesFlagsArrayName(int index)
495 {
496   return this->Internal->TK.getTimesFlagArray()[index].second.c_str();
497 }
498
499 int vtkMEDReader::GetTimesFlagsArrayStatus(const char *name)
500 {
501   if(!this->Internal)
502     return -1;
503   int pos(0);
504   std::istringstream iss(name); iss >> pos;
505   return (int)this->Internal->TK.getTimesFlagArray()[pos].first;
506 }
507
508 void vtkMEDReader::UpdateSIL(vtkInformation* request, vtkInformation *info)
509 {
510   if(!this->Internal)
511       return;
512   vtkMutableDirectedGraph *sil(vtkMutableDirectedGraph::New());
513
514   // This Should be more clever, TODO
515   std::string meshName(this->BuildSIL(sil));
516   if(meshName!=this->Internal->DftMeshName)
517     {
518       if(this->Internal->SIL)
519         this->Internal->SIL->Delete();
520       this->Internal->SIL=sil;
521       this->Internal->DftMeshName=meshName;
522     }
523   else
524     {
525       sil->Delete();
526     }
527 }
528
529 /*!
530  * The returned string is the name of the mesh activated which groups and families are in \a sil.
531  */
532 std::string vtkMEDReader::BuildSIL(vtkMutableDirectedGraph* sil)
533 {
534   if (!this->Internal)
535     return std::string();
536   sil->Initialize();
537   vtkSmartPointer<vtkVariantArray> childEdge(vtkSmartPointer<vtkVariantArray>::New());
538   childEdge->InsertNextValue(0);
539   vtkSmartPointer<vtkVariantArray> crossEdge(vtkSmartPointer<vtkVariantArray>::New());
540   crossEdge->InsertNextValue(1);
541   // CrossEdge is an edge linking hierarchies.
542   vtkUnsignedCharArray* crossEdgesArray=vtkUnsignedCharArray::New();
543   crossEdgesArray->SetName("CrossEdges");
544   sil->GetEdgeData()->AddArray(crossEdgesArray);
545   crossEdgesArray->Delete();
546   std::vector<std::string> names;
547   // Now build the hierarchy.
548   vtkIdType rootId=sil->AddVertex();
549   names.push_back("SIL");
550   // Add global fields root
551   vtkIdType fieldsRoot(sil->AddChild(rootId,childEdge));
552   names.push_back("FieldsStatusTree");
553   this->Internal->Tree.feedSIL(sil,fieldsRoot,childEdge,names);
554   vtkIdType meshesFamsGrpsRoot(sil->AddChild(rootId,childEdge));
555   names.push_back("MeshesFamsGrps");
556   std::string dftMeshName(this->Internal->Tree.feedSILForFamsAndGrps(sil,meshesFamsGrpsRoot,childEdge,names));
557   // This array is used to assign names to nodes.
558   vtkStringArray *namesArray(vtkStringArray::New());
559   namesArray->SetName("Names");
560   namesArray->SetNumberOfTuples(sil->GetNumberOfVertices());
561   sil->GetVertexData()->AddArray(namesArray);
562   namesArray->Delete();
563   std::vector<std::string>::const_iterator iter;
564   vtkIdType cc;
565   for(cc=0, iter=names.begin(); iter!=names.end(); ++iter, ++cc)
566     namesArray->SetValue(cc,(*iter).c_str());
567   return dftMeshName;
568 }
569
570 double vtkMEDReader::PublishTimeStepsIfNeeded(vtkInformation *outInfo, bool& isUpdated)
571 {
572   if(!this->Internal)
573     return 0.0;
574
575   int lev0(-1);
576   std::vector<double> tsteps;
577   if(!this->Internal->IsStdOrMode)
578     tsteps=this->Internal->Tree.getTimeSteps(lev0,this->Internal->TK);
579   else
580     { tsteps.resize(1); tsteps[0]=0.; }
581   isUpdated=false;
582   if(lev0!=this->Internal->LastLev0)
583     {
584       isUpdated=true;
585       double timeRange[2];
586       timeRange[0]=tsteps.front();
587       timeRange[1]=tsteps.back();
588       outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_STEPS(),&tsteps[0],tsteps.size());
589       outInfo->Set(vtkStreamingDemandDrivenPipeline::TIME_RANGE(),timeRange,2);
590       this->Internal->LastLev0=lev0;
591     }
592   return tsteps.front();
593 }
594
595 void vtkMEDReader::FillMultiBlockDataSetInstance(vtkMultiBlockDataSet *output, double reqTS)
596 {
597   if( !this->Internal )
598     return;
599   std::string meshName;
600   vtkDataSet *ret(this->Internal->Tree.buildVTKInstance(this->Internal->IsStdOrMode,reqTS,meshName,this->Internal->TK));
601   if(this->Internal->GenerateVect)
602     {
603       vtkGenerateVectors::Operate(ret->GetPointData());
604       vtkGenerateVectors::Operate(ret->GetCellData());
605       vtkGenerateVectors::Operate(ret->GetFieldData());
606     }
607   output->SetBlock(0,ret);
608   ret->Delete();
609 }
610
611 void vtkMEDReader::PrintSelf(ostream& os, vtkIndent indent)
612 {
613   this->Superclass::PrintSelf(os, indent);
614 }