Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-04-10 08:06:20

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 #include "DicomFilePlan.hh"
0027 
0028 #include "DicomBeam.hh"
0029 #include "DicomBeamBlock.hh"
0030 #include "DicomBeamCompensator.hh"
0031 #include "DicomBeamControlPoint.hh"
0032 #include "DicomBeamDevicePos.hh"
0033 #include "DicomBeamDeviceRef.hh"
0034 #include "DicomBeamWedge.hh"
0035 #include "dcmtk/dcmdata/dcdeftag.h"
0036 #include "dcmtk/dcmdata/dcfilefo.h"
0037 #include "dcmtk/dcmrt/drtplan.h"
0038 #include "dcmtk/dcmrt/seq/drtbl2.h"  // for BlockSequence
0039 #include "dcmtk/dcmrt/seq/drtbldps.h"  // for BeamLimitingDevicePositionSequence
0040 #include "dcmtk/dcmrt/seq/drtblds1.h"  // for BeamLimitingDeviceSequence
0041 #include "dcmtk/dcmrt/seq/drtbs.h"  // for BeamSequence
0042 #include "dcmtk/dcmrt/seq/drtcos.h"  // for CompensatorSequence
0043 #include "dcmtk/dcmrt/seq/drtcps.h"  // for ControlPointSequence
0044 #include "dcmtk/dcmrt/seq/drtfgs.h"  // DRTFractionGroupSequence
0045 #include "dcmtk/dcmrt/seq/drtrbs8.h"  // DRTReferencedBeamSequenceInRTFractionSchemeModule
0046 #include "dcmtk/dcmrt/seq/drtws.h"  // for WedgeSequence
0047 
0048 #include "G4ThreeVector.hh"
0049 
0050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0051 DicomFilePlan::DicomFilePlan(DcmDataset* dset) : DicomVFile(dset) {}
0052 
0053 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0054 void DicomFilePlan::ReadData()
0055 {
0056   DRTPlanIOD rtplan;
0057   OFCondition result = rtplan.read(*theDataset);
0058   if (!result.good()) {
0059     G4Exception("DicomFilePlan::ReadData", "DFS001", FatalException, result.text());
0060   }
0061   OFString fstr;
0062   Sint32 fint;
0063   Float64 ffloat;
0064   OFVector<Float64> fvfloat;
0065 
0066   DRTFractionGroupSequence frgSeq = rtplan.getFractionGroupSequence();
0067   if (frgSeq.isEmpty()) {
0068     G4Exception("DicomFilePlan::ReadData", "DFS002", JustWarning,
0069                 "DRTFractionGroupSequence is empty");
0070   }
0071   G4cout << "@@@@@ NUMBER OF FractionSequences " << frgSeq.getNumberOfItems() << G4endl;
0072   frgSeq.gotoFirstItem();
0073   for (size_t i1 = 0; i1 < frgSeq.getNumberOfItems(); i1++) {
0074     DRTFractionGroupSequence::Item& rfgItem = frgSeq.getCurrentItem();
0075     rfgItem.getBeamDoseMeaning(fstr);
0076     G4cout << " " << i1 << " BeamDoseMeaning " << fstr << G4endl;
0077     rfgItem.getFractionGroupDescription(fstr);
0078     G4cout << " " << i1 << " FractionGroupDescription " << fstr << G4endl;
0079     rfgItem.getFractionGroupNumber(fint);
0080     G4cout << " " << i1 << " FractionGroupNumber " << fint << G4endl;
0081     rfgItem.getFractionPattern(fstr);
0082     G4cout << " " << i1 << " FractionPattern " << fstr << G4endl;
0083     rfgItem.getNumberOfBeams(fint);
0084     G4cout << " " << i1 << " NumberOfBeams " << fint << G4endl;
0085     theNumberOfBeams = fint;
0086     rfgItem.getNumberOfBrachyApplicationSetups(fint);
0087     G4cout << " " << i1 << " NumberOfBrachyApplicationSetups " << fint << G4endl;
0088     CheckData0(" NumberOfBrachyApplicationSetups ", fint);
0089     rfgItem.getNumberOfFractionPatternDigitsPerDay(fint);
0090     G4cout << " " << i1 << " NumberOfFractionPatternDigitsPerDay " << fint << G4endl;
0091     rfgItem.getNumberOfFractionsPlanned(fint);
0092     G4cout << " " << i1 << " NumberOfFractionsPlanned " << fint << G4endl;
0093     rfgItem.getRepeatFractionCycleLength(fint);
0094     G4cout << " " << i1 << " RepeatFractionCycleLength " << fint << G4endl;
0095     DRTReferencedBeamSequenceInRTFractionSchemeModule refBeamSeq =
0096       rfgItem.getReferencedBeamSequence();
0097     G4cout << " @@@ NUMBER OF ReferencedBeamSequences " << refBeamSeq.getNumberOfItems() << G4endl;
0098 
0099     refBeamSeq.gotoFirstItem();
0100     for (size_t i2 = 0; i2 < refBeamSeq.getNumberOfItems(); i2++) {
0101       DicomBeam* db = new DicomBeam();
0102       theBeams.push_back(db);
0103       DRTReferencedBeamSequenceInRTFractionSchemeModule::Item& rbsItem =
0104         refBeamSeq.getCurrentItem();
0105       rbsItem.getBeamDeliveryDurationLimit(ffloat);
0106       G4cout << "  " << i2 << " BeamDeliveryDurationLimit " << ffloat << G4endl;
0107       rbsItem.getBeamDose(ffloat);
0108       G4cout << "  " << i2 << " BeamDose " << ffloat << G4endl;  // dose at dose point
0109       rbsItem.getBeamDoseSpecificationPoint(fvfloat);
0110       G4cout << "  " << i2 << " BeamDoseSpecificationPoint (" << fvfloat[0] << "," << fvfloat[1]
0111              << "," << fvfloat[2] << ")" << G4endl;
0112       db->SetDoseSpecificationPoint(G4ThreeVector(fvfloat[0], fvfloat[1], fvfloat[2]));
0113       rbsItem.getBeamMeterset(ffloat);
0114       G4cout << "  " << i2 << " BeamMeterset " << ffloat << G4endl;
0115       db->SetMeterset(ffloat);
0116       rbsItem.getReferencedBeamNumber(fint);
0117       G4cout << "  " << i2 << " ReferencedBeamNumber " << fint << G4endl;
0118 
0119       refBeamSeq.gotoNextItem();
0120     }
0121 
0122     frgSeq.gotoNextItem();
0123   }
0124 
0125   DRTBeamSequence beamSeq = rtplan.getBeamSequence();
0126   if (beamSeq.isEmpty()) {
0127     G4Exception("DicomFilePlan::ReadData", "DFS002", JustWarning, "DRTBeamSequence is empty");
0128   }
0129   G4cout << "@@@@@ NUMBER OF BeamSequences " << beamSeq.getNumberOfItems() << G4endl;
0130   beamSeq.gotoFirstItem();
0131   for (size_t i1 = 0; i1 < beamSeq.getNumberOfItems(); i1++) {
0132     DicomBeam* db = theBeams[i1];
0133     DRTBeamSequence::Item& beamItem = beamSeq.getCurrentItem();
0134 
0135     beamItem.getManufacturer(fstr);
0136     G4cout << " " << i1 << " Manufacturer " << fstr << G4endl;
0137     beamItem.getManufacturerModelName(fstr);
0138     G4cout << " " << i1 << " ManufacturerModelName " << fstr << G4endl;
0139     beamItem.getTreatmentMachineName(fstr);
0140     G4cout << " " << i1 << " TreatmentMachineName " << fstr << G4endl;
0141     beamItem.getPrimaryDosimeterUnit(fstr);
0142     G4cout << " " << i1 << " PrimaryDosimeterUnit " << fstr << G4endl;
0143     beamItem.getSourceAxisDistance(ffloat);
0144     G4cout << " " << i1 << " SourceAxisDistance " << ffloat << G4endl;
0145     db->SetSourceAxisDistance(ffloat);
0146 
0147     DRTBeamLimitingDeviceSequenceInRTBeamsModule beamLDS = beamItem.getBeamLimitingDeviceSequence();
0148     G4cout << " @@@ NUMBER OF BeamLimitingDeviceSequence " << beamLDS.getNumberOfItems() << G4endl;
0149     beamLDS.gotoFirstItem();
0150     for (size_t i2 = 0; i2 < beamLDS.getNumberOfItems(); i2++) {
0151       DRTBeamLimitingDeviceSequenceInRTBeamsModule::Item bldsItem = beamLDS.getCurrentItem();
0152       DicomBeamDeviceRef* dbd = new DicomBeamDeviceRef(bldsItem);
0153       db->AddDevice(dbd);
0154 
0155       beamLDS.gotoNextItem();
0156     }
0157 
0158     beamItem.getBeamNumber(fint);
0159     G4cout << " " << i1 << " BeamNumber " << fint << G4endl;
0160     db->SetNumber(fint);
0161     beamItem.getBeamName(fstr);
0162     G4cout << " " << i1 << " BeamName " << fstr << G4endl;
0163     beamItem.getBeamDescription(fstr);
0164     G4cout << " " << i1 << " BeamDescription " << fstr << G4endl;
0165     beamItem.getBeamType(fstr);
0166     G4cout << " " << i1 << " BeamType " << fstr << G4endl;
0167     beamItem.getRadiationType(fstr);
0168     G4cout << " " << i1 << " RadiationType " << fstr << G4endl;
0169     db->SetRadiationType(fstr);
0170     beamItem.getTreatmentDeliveryType(fstr);
0171     G4cout << " " << i1 << " TreatmentDeliveryType " << fstr << G4endl;
0172     beamItem.getNumberOfWedges(fint);
0173     G4cout << " " << i1 << " NumberOfWedges " << fint << G4endl;
0174     DRTWedgeSequence beamWedge = beamItem.getWedgeSequence();
0175     beamWedge.gotoFirstItem();
0176     for (size_t i2 = 0; i2 < beamWedge.getNumberOfItems(); i2++) {
0177       DRTWedgeSequence::Item bwedItem = beamWedge.getCurrentItem();
0178       DicomBeamWedge* dbwed = new DicomBeamWedge(bwedItem);
0179       db->AddWedge(dbwed);
0180       beamWedge.gotoNextItem();
0181     }
0182 
0183     beamItem.getNumberOfCompensators(fint);
0184     G4cout << " " << i1 << " NumberOfCompensators " << fint << G4endl;
0185     DRTCompensatorSequence beamCompens = beamItem.getCompensatorSequence();
0186     beamCompens.gotoFirstItem();
0187     for (size_t i2 = 0; i2 < beamCompens.getNumberOfItems(); i2++) {
0188       DRTCompensatorSequence::Item bcompItem = beamCompens.getCurrentItem();
0189       DicomBeamCompensator* dbcomp = new DicomBeamCompensator(bcompItem);
0190       db->AddCompensator(dbcomp);
0191       beamCompens.gotoNextItem();
0192     }
0193 
0194     beamItem.getNumberOfBoli(fint);
0195     G4cout << " " << i1 << " NumberOfBoli " << fint << G4endl;
0196     // Bolus has no relevant info (see drtrbos1.h)
0197 
0198     beamItem.getNumberOfBlocks(fint);
0199     G4cout << " " << i1 << " NumberOfBlocks " << fint << G4endl;
0200     DRTBlockSequenceInRTBeamsModule beamBlock = beamItem.getBlockSequence();
0201     beamBlock.gotoFirstItem();
0202     for (size_t i2 = 0; i2 < beamBlock.getNumberOfItems(); i2++) {
0203       DRTBlockSequenceInRTBeamsModule::Item bblItem = beamBlock.getCurrentItem();
0204       DicomBeamBlock* dbbl = new DicomBeamBlock(bblItem);
0205       db->AddBlock(dbbl);
0206       beamBlock.gotoNextItem();
0207     }
0208 
0209     beamItem.getFinalCumulativeMetersetWeight(fstr);
0210     G4cout << " " << i1 << " FinalCumulativeMetersetWeight " << fstr << G4endl;
0211     beamItem.getDeviceSerialNumber(fstr);
0212     G4cout << " " << i1 << " DeviceSerialNumber " << fstr << G4endl;
0213     beamItem.getHighDoseTechniqueType(fstr);
0214     G4cout << " " << i1 << " HighDoseTechniqueType " << fstr << G4endl;
0215     beamItem.getInstitutionAddress(fstr);
0216     G4cout << " " << i1 << " InstitutionAddress " << fstr << G4endl;
0217     beamItem.getInstitutionName(fstr);
0218     G4cout << " " << i1 << " InstitutionName " << fstr << G4endl;
0219     beamItem.getInstitutionalDepartmentName(fstr);
0220     G4cout << " " << i1 << " InstitutionalDepartmentName " << fstr << G4endl;
0221     beamItem.getReferencedPatientSetupNumber(fint);
0222     G4cout << " " << i1 << " ReferencedPatientSetupNumber " << fint << G4endl;
0223     beamItem.getReferencedToleranceTableNumber(fint);
0224     G4cout << " " << i1 << " ReferencedToleranceTableNumber " << fint << G4endl;
0225     beamItem.getTotalBlockTrayFactor(ffloat);
0226     G4cout << " " << i1 << " TotalBlockTrayFactor " << ffloat << G4endl;
0227     beamItem.getTotalCompensatorTrayFactor(ffloat);
0228     G4cout << " " << i1 << " TotalCompensatorTrayFactor " << ffloat << G4endl;
0229 
0230     beamItem.getNumberOfControlPoints(fint);
0231     DRTControlPointSequence controlPSeq = beamItem.getControlPointSequence();
0232     G4cout << " @@@ NUMBER OF ControlPointSequences " << controlPSeq.getNumberOfItems() << " = "
0233            << fint << G4endl;
0234     controlPSeq.gotoFirstItem();
0235     DicomBeamControlPoint* dbcp0 = 0;
0236     // Only first ControlPoint has some info if it does not change
0237     for (size_t i2 = 0; i2 < controlPSeq.getNumberOfItems(); i2++) {
0238       DRTControlPointSequence::Item& cpItem = controlPSeq.getCurrentItem();
0239       if (db->GetNControlPoints() != 0) dbcp0 = db->GetControlPoint(0);
0240       DicomBeamControlPoint* dbcp = new DicomBeamControlPoint(cpItem, dbcp0);
0241       db->AddControlPoint(dbcp);
0242       controlPSeq.gotoNextItem();
0243     }
0244 
0245     beamSeq.gotoNextItem();
0246   }
0247 
0248   //(300a,0180) SQ (Sequence with explicit length #=1)      #  30, 1 PatientSetupSequence
0249   //(300c,0060) SQ (Sequence with explicit length #=1)      # 112, 1 ReferencedStructureSetSequence
0250   //(300e,0002) CS [UNAPPROVED]                             #  10, 1 ApprovalStatus
0251   //(7fe0,0010) OW (no value available)                     #   0, 1 PixelData
0252 }
0253 
0254 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0255 void DicomFilePlan::CheckData0(OFString title, Sint32 val)
0256 {
0257   if (val != 0) {
0258     G4Exception("DicomFilePlan::CheckData", "DFP003", FatalException,
0259                 (title + " exists, and code is not ready ").c_str());
0260   }
0261 }
0262 
0263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0264 void DicomFilePlan::SetControlPointMetersets()
0265 {
0266   for (size_t ii = 0; ii < theBeams.size(); ii++) {
0267     theBeams[ii]->SetControlPointMetersets();
0268   }
0269 }
0270 
0271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
0272 void DicomFilePlan::DumpToFile()
0273 {
0274   for (size_t ii = 0; ii < theBeams.size(); ii++) {
0275     theBeams[ii]->DumpToFile();
0276   }
0277 }