Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/examples/extended/medical/DICOM/dicomReader/src/DicomFilePlan.cc

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  1 //
  2 // ********************************************************************
  3 // * License and Disclaimer                                           *
  4 // *                                                                  *
  5 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
  6 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
  7 // * conditions of the Geant4 Software License,  included in the file *
  8 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
  9 // * include a list of copyright holders.                             *
 10 // *                                                                  *
 11 // * Neither the authors of this software system, nor their employing *
 12 // * institutes,nor the agencies providing financial support for this *
 13 // * work  make  any representation or  warranty, express or implied, *
 14 // * regarding  this  software system or assume any liability for its *
 15 // * use.  Please see the license in the file  LICENSE  and URL above *
 16 // * for the full disclaimer and the limitation of liability.         *
 17 // *                                                                  *
 18 // * This  code  implementation is the result of  the  scientific and *
 19 // * technical work of the GEANT4 collaboration.                      *
 20 // * By using,  copying,  modifying or  distributing the software (or *
 21 // * any work based  on the software)  you  agree  to acknowledge its *
 22 // * use  in  resulting  scientific  publications,  and indicate your *
 23 // * acceptance of all terms of the Geant4 Software license.          *
 24 // ********************************************************************
 25 //
 26 #include "DicomFilePlan.hh"
 27 
 28 #include "DicomBeam.hh"
 29 #include "DicomBeamBlock.hh"
 30 #include "DicomBeamCompensator.hh"
 31 #include "DicomBeamControlPoint.hh"
 32 #include "DicomBeamDevicePos.hh"
 33 #include "DicomBeamDeviceRef.hh"
 34 #include "DicomBeamWedge.hh"
 35 #include "dcmtk/dcmdata/dcdeftag.h"
 36 #include "dcmtk/dcmdata/dcfilefo.h"
 37 #include "dcmtk/dcmrt/drtplan.h"
 38 #include "dcmtk/dcmrt/seq/drtbl2.h"  // for BlockSequence
 39 #include "dcmtk/dcmrt/seq/drtbldps.h"  // for BeamLimitingDevicePositionSequence
 40 #include "dcmtk/dcmrt/seq/drtblds1.h"  // for BeamLimitingDeviceSequence
 41 #include "dcmtk/dcmrt/seq/drtbs.h"  // for BeamSequence
 42 #include "dcmtk/dcmrt/seq/drtcos.h"  // for CompensatorSequence
 43 #include "dcmtk/dcmrt/seq/drtcps.h"  // for ControlPointSequence
 44 #include "dcmtk/dcmrt/seq/drtfgs.h"  // DRTFractionGroupSequence
 45 #include "dcmtk/dcmrt/seq/drtrbs8.h"  // DRTReferencedBeamSequenceInRTFractionSchemeModule
 46 #include "dcmtk/dcmrt/seq/drtws.h"  // for WedgeSequence
 47 
 48 #include "G4ThreeVector.hh"
 49 
 50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 51 DicomFilePlan::DicomFilePlan(DcmDataset* dset) : DicomVFile(dset) {}
 52 
 53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
 54 void DicomFilePlan::ReadData()
 55 {
 56   DRTPlanIOD rtplan;
 57   OFCondition result = rtplan.read(*theDataset);
 58   if (!result.good()) {
 59     G4Exception("DicomFilePlan::ReadData", "DFS001", FatalException, result.text());
 60   }
 61   OFString fstr;
 62   Sint32 fint;
 63   Float64 ffloat;
 64   OFVector<Float64> fvfloat;
 65 
 66   DRTFractionGroupSequence frgSeq = rtplan.getFractionGroupSequence();
 67   if (frgSeq.isEmpty()) {
 68     G4Exception("DicomFilePlan::ReadData", "DFS002", JustWarning,
 69                 "DRTFractionGroupSequence is empty");
 70   }
 71   G4cout << "@@@@@ NUMBER OF FractionSequences " << frgSeq.getNumberOfItems() << G4endl;
 72   frgSeq.gotoFirstItem();
 73   for (size_t i1 = 0; i1 < frgSeq.getNumberOfItems(); i1++) {
 74     DRTFractionGroupSequence::Item& rfgItem = frgSeq.getCurrentItem();
 75     rfgItem.getBeamDoseMeaning(fstr);
 76     G4cout << " " << i1 << " BeamDoseMeaning " << fstr << G4endl;
 77     rfgItem.getFractionGroupDescription(fstr);
 78     G4cout << " " << i1 << " FractionGroupDescription " << fstr << G4endl;
 79     rfgItem.getFractionGroupNumber(fint);
 80     G4cout << " " << i1 << " FractionGroupNumber " << fint << G4endl;
 81     rfgItem.getFractionPattern(fstr);
 82     G4cout << " " << i1 << " FractionPattern " << fstr << G4endl;
 83     rfgItem.getNumberOfBeams(fint);
 84     G4cout << " " << i1 << " NumberOfBeams " << fint << G4endl;
 85     theNumberOfBeams = fint;
 86     rfgItem.getNumberOfBrachyApplicationSetups(fint);
 87     G4cout << " " << i1 << " NumberOfBrachyApplicationSetups " << fint << G4endl;
 88     CheckData0(" NumberOfBrachyApplicationSetups ", fint);
 89     rfgItem.getNumberOfFractionPatternDigitsPerDay(fint);
 90     G4cout << " " << i1 << " NumberOfFractionPatternDigitsPerDay " << fint << G4endl;
 91     rfgItem.getNumberOfFractionsPlanned(fint);
 92     G4cout << " " << i1 << " NumberOfFractionsPlanned " << fint << G4endl;
 93     rfgItem.getRepeatFractionCycleLength(fint);
 94     G4cout << " " << i1 << " RepeatFractionCycleLength " << fint << G4endl;
 95     DRTReferencedBeamSequenceInRTFractionSchemeModule refBeamSeq =
 96       rfgItem.getReferencedBeamSequence();
 97     G4cout << " @@@ NUMBER OF ReferencedBeamSequences " << refBeamSeq.getNumberOfItems() << G4endl;
 98 
 99     refBeamSeq.gotoFirstItem();
100     for (size_t i2 = 0; i2 < refBeamSeq.getNumberOfItems(); i2++) {
101       DicomBeam* db = new DicomBeam();
102       theBeams.push_back(db);
103       DRTReferencedBeamSequenceInRTFractionSchemeModule::Item& rbsItem =
104         refBeamSeq.getCurrentItem();
105       rbsItem.getBeamDeliveryDurationLimit(ffloat);
106       G4cout << "  " << i2 << " BeamDeliveryDurationLimit " << ffloat << G4endl;
107       rbsItem.getBeamDose(ffloat);
108       G4cout << "  " << i2 << " BeamDose " << ffloat << G4endl;  // dose at dose point
109       rbsItem.getBeamDoseSpecificationPoint(fvfloat);
110       G4cout << "  " << i2 << " BeamDoseSpecificationPoint (" << fvfloat[0] << "," << fvfloat[1]
111              << "," << fvfloat[2] << ")" << G4endl;
112       db->SetDoseSpecificationPoint(G4ThreeVector(fvfloat[0], fvfloat[1], fvfloat[2]));
113       rbsItem.getBeamMeterset(ffloat);
114       G4cout << "  " << i2 << " BeamMeterset " << ffloat << G4endl;
115       db->SetMeterset(ffloat);
116       rbsItem.getReferencedBeamNumber(fint);
117       G4cout << "  " << i2 << " ReferencedBeamNumber " << fint << G4endl;
118 
119       refBeamSeq.gotoNextItem();
120     }
121 
122     frgSeq.gotoNextItem();
123   }
124 
125   DRTBeamSequence beamSeq = rtplan.getBeamSequence();
126   if (beamSeq.isEmpty()) {
127     G4Exception("DicomFilePlan::ReadData", "DFS002", JustWarning, "DRTBeamSequence is empty");
128   }
129   G4cout << "@@@@@ NUMBER OF BeamSequences " << beamSeq.getNumberOfItems() << G4endl;
130   beamSeq.gotoFirstItem();
131   for (size_t i1 = 0; i1 < beamSeq.getNumberOfItems(); i1++) {
132     DicomBeam* db = theBeams[i1];
133     DRTBeamSequence::Item& beamItem = beamSeq.getCurrentItem();
134 
135     beamItem.getManufacturer(fstr);
136     G4cout << " " << i1 << " Manufacturer " << fstr << G4endl;
137     beamItem.getManufacturerModelName(fstr);
138     G4cout << " " << i1 << " ManufacturerModelName " << fstr << G4endl;
139     beamItem.getTreatmentMachineName(fstr);
140     G4cout << " " << i1 << " TreatmentMachineName " << fstr << G4endl;
141     beamItem.getPrimaryDosimeterUnit(fstr);
142     G4cout << " " << i1 << " PrimaryDosimeterUnit " << fstr << G4endl;
143     beamItem.getSourceAxisDistance(ffloat);
144     G4cout << " " << i1 << " SourceAxisDistance " << ffloat << G4endl;
145     db->SetSourceAxisDistance(ffloat);
146 
147     DRTBeamLimitingDeviceSequenceInRTBeamsModule beamLDS = beamItem.getBeamLimitingDeviceSequence();
148     G4cout << " @@@ NUMBER OF BeamLimitingDeviceSequence " << beamLDS.getNumberOfItems() << G4endl;
149     beamLDS.gotoFirstItem();
150     for (size_t i2 = 0; i2 < beamLDS.getNumberOfItems(); i2++) {
151       DRTBeamLimitingDeviceSequenceInRTBeamsModule::Item bldsItem = beamLDS.getCurrentItem();
152       DicomBeamDeviceRef* dbd = new DicomBeamDeviceRef(bldsItem);
153       db->AddDevice(dbd);
154 
155       beamLDS.gotoNextItem();
156     }
157 
158     beamItem.getBeamNumber(fint);
159     G4cout << " " << i1 << " BeamNumber " << fint << G4endl;
160     db->SetNumber(fint);
161     beamItem.getBeamName(fstr);
162     G4cout << " " << i1 << " BeamName " << fstr << G4endl;
163     beamItem.getBeamDescription(fstr);
164     G4cout << " " << i1 << " BeamDescription " << fstr << G4endl;
165     beamItem.getBeamType(fstr);
166     G4cout << " " << i1 << " BeamType " << fstr << G4endl;
167     beamItem.getRadiationType(fstr);
168     G4cout << " " << i1 << " RadiationType " << fstr << G4endl;
169     db->SetRadiationType(fstr);
170     beamItem.getTreatmentDeliveryType(fstr);
171     G4cout << " " << i1 << " TreatmentDeliveryType " << fstr << G4endl;
172     beamItem.getNumberOfWedges(fint);
173     G4cout << " " << i1 << " NumberOfWedges " << fint << G4endl;
174     DRTWedgeSequence beamWedge = beamItem.getWedgeSequence();
175     beamWedge.gotoFirstItem();
176     for (size_t i2 = 0; i2 < beamWedge.getNumberOfItems(); i2++) {
177       DRTWedgeSequence::Item bwedItem = beamWedge.getCurrentItem();
178       DicomBeamWedge* dbwed = new DicomBeamWedge(bwedItem);
179       db->AddWedge(dbwed);
180       beamWedge.gotoNextItem();
181     }
182 
183     beamItem.getNumberOfCompensators(fint);
184     G4cout << " " << i1 << " NumberOfCompensators " << fint << G4endl;
185     DRTCompensatorSequence beamCompens = beamItem.getCompensatorSequence();
186     beamCompens.gotoFirstItem();
187     for (size_t i2 = 0; i2 < beamCompens.getNumberOfItems(); i2++) {
188       DRTCompensatorSequence::Item bcompItem = beamCompens.getCurrentItem();
189       DicomBeamCompensator* dbcomp = new DicomBeamCompensator(bcompItem);
190       db->AddCompensator(dbcomp);
191       beamCompens.gotoNextItem();
192     }
193 
194     beamItem.getNumberOfBoli(fint);
195     G4cout << " " << i1 << " NumberOfBoli " << fint << G4endl;
196     // Bolus has no relevant info (see drtrbos1.h)
197 
198     beamItem.getNumberOfBlocks(fint);
199     G4cout << " " << i1 << " NumberOfBlocks " << fint << G4endl;
200     DRTBlockSequenceInRTBeamsModule beamBlock = beamItem.getBlockSequence();
201     beamBlock.gotoFirstItem();
202     for (size_t i2 = 0; i2 < beamBlock.getNumberOfItems(); i2++) {
203       DRTBlockSequenceInRTBeamsModule::Item bblItem = beamBlock.getCurrentItem();
204       DicomBeamBlock* dbbl = new DicomBeamBlock(bblItem);
205       db->AddBlock(dbbl);
206       beamBlock.gotoNextItem();
207     }
208 
209     beamItem.getFinalCumulativeMetersetWeight(fstr);
210     G4cout << " " << i1 << " FinalCumulativeMetersetWeight " << fstr << G4endl;
211     beamItem.getDeviceSerialNumber(fstr);
212     G4cout << " " << i1 << " DeviceSerialNumber " << fstr << G4endl;
213     beamItem.getHighDoseTechniqueType(fstr);
214     G4cout << " " << i1 << " HighDoseTechniqueType " << fstr << G4endl;
215     beamItem.getInstitutionAddress(fstr);
216     G4cout << " " << i1 << " InstitutionAddress " << fstr << G4endl;
217     beamItem.getInstitutionName(fstr);
218     G4cout << " " << i1 << " InstitutionName " << fstr << G4endl;
219     beamItem.getInstitutionalDepartmentName(fstr);
220     G4cout << " " << i1 << " InstitutionalDepartmentName " << fstr << G4endl;
221     beamItem.getReferencedPatientSetupNumber(fint);
222     G4cout << " " << i1 << " ReferencedPatientSetupNumber " << fint << G4endl;
223     beamItem.getReferencedToleranceTableNumber(fint);
224     G4cout << " " << i1 << " ReferencedToleranceTableNumber " << fint << G4endl;
225     beamItem.getTotalBlockTrayFactor(ffloat);
226     G4cout << " " << i1 << " TotalBlockTrayFactor " << ffloat << G4endl;
227     beamItem.getTotalCompensatorTrayFactor(ffloat);
228     G4cout << " " << i1 << " TotalCompensatorTrayFactor " << ffloat << G4endl;
229 
230     beamItem.getNumberOfControlPoints(fint);
231     DRTControlPointSequence controlPSeq = beamItem.getControlPointSequence();
232     G4cout << " @@@ NUMBER OF ControlPointSequences " << controlPSeq.getNumberOfItems() << " = "
233            << fint << G4endl;
234     controlPSeq.gotoFirstItem();
235     DicomBeamControlPoint* dbcp0 = 0;
236     // Only first ControlPoint has some info if it does not change
237     for (size_t i2 = 0; i2 < controlPSeq.getNumberOfItems(); i2++) {
238       DRTControlPointSequence::Item& cpItem = controlPSeq.getCurrentItem();
239       if (db->GetNControlPoints() != 0) dbcp0 = db->GetControlPoint(0);
240       DicomBeamControlPoint* dbcp = new DicomBeamControlPoint(cpItem, dbcp0);
241       db->AddControlPoint(dbcp);
242       controlPSeq.gotoNextItem();
243     }
244 
245     beamSeq.gotoNextItem();
246   }
247 
248   //(300a,0180) SQ (Sequence with explicit length #=1)      #  30, 1 PatientSetupSequence
249   //(300c,0060) SQ (Sequence with explicit length #=1)      # 112, 1 ReferencedStructureSetSequence
250   //(300e,0002) CS [UNAPPROVED]                             #  10, 1 ApprovalStatus
251   //(7fe0,0010) OW (no value available)                     #   0, 1 PixelData
252 }
253 
254 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
255 void DicomFilePlan::CheckData0(OFString title, Sint32 val)
256 {
257   if (val != 0) {
258     G4Exception("DicomFilePlan::CheckData", "DFP003", FatalException,
259                 (title + " exists, and code is not ready ").c_str());
260   }
261 }
262 
263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
264 void DicomFilePlan::SetControlPointMetersets()
265 {
266   for (size_t ii = 0; ii < theBeams.size(); ii++) {
267     theBeams[ii]->SetControlPointMetersets();
268   }
269 }
270 
271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
272 void DicomFilePlan::DumpToFile()
273 {
274   for (size_t ii = 0; ii < theBeams.size(); ii++) {
275     theBeams[ii]->DumpToFile();
276   }
277 }
278