matRad_daoVec2ApertureInfo

Purpose ^

matRad function to translate vector representation into struct

Synopsis ^

function updatedInfo = matRad_daoVec2ApertureInfo(apertureInfo,apertureInfoVect)

Description ^

 matRad function to translate vector representation into struct
 The vector representation of the aperture shape and weight are translated 
 into an aperture info struct. At the same time, the updated bixel weight
 vector w is computed and a vector listing the correspondence between leaf 
 tips and bixel indices for gradient calculation

 call
   updatedInfo = matRad_daoVec2ApertureInfo(apertureInfo,apertureInfoVect)

 input
   apertureInfo:     aperture shape info struct
   apertureInfoVect: aperture weights and shapes parameterized as vector

 output
   updatedInfo: updated aperture shape info struct according to apertureInfoVect

 References
   

 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%

Cross-reference information ^

This function calls: This function is called by:

Source code ^

0001 function updatedInfo = matRad_daoVec2ApertureInfo(apertureInfo,apertureInfoVect)
0002 % matRad function to translate vector representation into struct
0003 % The vector representation of the aperture shape and weight are translated
0004 % into an aperture info struct. At the same time, the updated bixel weight
0005 % vector w is computed and a vector listing the correspondence between leaf
0006 % tips and bixel indices for gradient calculation
0007 %
0008 % call
0009 %   updatedInfo = matRad_daoVec2ApertureInfo(apertureInfo,apertureInfoVect)
0010 %
0011 % input
0012 %   apertureInfo:     aperture shape info struct
0013 %   apertureInfoVect: aperture weights and shapes parameterized as vector
0014 %
0015 % output
0016 %   updatedInfo: updated aperture shape info struct according to apertureInfoVect
0017 %
0018 % References
0019 %
0020 %
0021 % %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0022 
0023 % %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0024 %
0025 % Copyright 2015 the matRad development team.
0026 %
0027 % This file is part of the matRad project. It is subject to the license
0028 % terms in the LICENSE file found in the top-level directory of this
0029 % distribution and at https://github.com/e0404/matRad/LICENSES.txt. No part
0030 % of the matRad project, including this file, may be copied, modified,
0031 % propagated, or distributed except according to the terms contained in the
0032 % LICENSE file.
0033 %
0034 % %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
0035 
0036 % function to update the apertureInfo struct after the each iteraton of the
0037 % optimization
0038 
0039 w = zeros(apertureInfo.totalNumOfBixels,1);
0040 
0041 % initializing variables
0042 updatedInfo = apertureInfo;
0043 
0044 updatedInfo.apertureVector = apertureInfoVect;
0045 
0046 shapeInd = 1;
0047 
0048 indVect = NaN*ones(apertureInfo.totalNumOfShapes + apertureInfo.totalNumOfLeafPairs,1);
0049 
0050 % helper function to cope with numerical instabilities through rounding
0051 round2 = @(a,b) round(a*10^b)/10^b;
0052 
0053 
0054 %% update the shapeMaps
0055 % here the new colimator positions are used to create new shapeMaps that
0056 % now include decimal values instead of binary
0057 
0058 % loop over all beams
0059 for i = 1:numel(updatedInfo.beam)
0060     
0061     %posOfRightCornerPixel = apertureInfo.beam(i).posOfCornerBixel(1) + (size(apertureInfo.beam(i).bixelIndMap,2)-1)*apertureInfo.bixelWidth;
0062 
0063     % pre compute left and right bixel edges
0064     edges_l = updatedInfo.beam(i).posOfCornerBixel(1)...
0065                  + ([1:size(apertureInfo.beam(i).bixelIndMap,2)]-1-1/2)*updatedInfo.bixelWidth;
0066     edges_r = updatedInfo.beam(i).posOfCornerBixel(1)...
0067                  + ([1:size(apertureInfo.beam(i).bixelIndMap,2)]-1+1/2)*updatedInfo.bixelWidth;
0068     
0069     % loop over all shapes
0070     for j = 1:updatedInfo.beam(i).numOfShapes
0071         
0072         % update the shape weight
0073         updatedInfo.beam(i).shape(j).weight = apertureInfoVect(shapeInd);
0074         
0075         % get dimensions of 2d matrices that store shape/bixel information
0076         n = apertureInfo.beam(i).numOfActiveLeafPairs;
0077         
0078         % extract left and right leaf positions from shape vector
0079         vectorIx     = updatedInfo.beam(i).shape(j).vectorOffset + ([1:n]-1);
0080         leftLeafPos  = apertureInfoVect(vectorIx);
0081         rightLeafPos = apertureInfoVect(vectorIx+apertureInfo.totalNumOfLeafPairs);
0082         
0083         % update information in shape structure
0084         updatedInfo.beam(i).shape(j).leftLeafPos  = leftLeafPos;
0085         updatedInfo.beam(i).shape(j).rightLeafPos = rightLeafPos;
0086         
0087         % rounding for numerical stability
0088         leftLeafPos  = round2(leftLeafPos,6);
0089         rightLeafPos = round2(rightLeafPos,6);
0090         
0091         %
0092         xPosIndLeftLeaf  = round((leftLeafPos - apertureInfo.beam(i).posOfCornerBixel(1))/apertureInfo.bixelWidth + 1);
0093         xPosIndRightLeaf = round((rightLeafPos - apertureInfo.beam(i).posOfCornerBixel(1))/apertureInfo.bixelWidth + 1);
0094         
0095         % check limits because of rounding off issues at maximum, i.e.,
0096         % enforce round(X.5) -> X
0097         xPosIndLeftLeaf(leftLeafPos == apertureInfo.beam(i).lim_r) = ...
0098             .5 + (leftLeafPos(leftLeafPos == apertureInfo.beam(i).lim_r) ...
0099             - apertureInfo.beam(i).posOfCornerBixel(1))/apertureInfo.bixelWidth;
0100         xPosIndRightLeaf(rightLeafPos == apertureInfo.beam(i).lim_r) = ...
0101             .5 + (rightLeafPos(rightLeafPos == apertureInfo.beam(i).lim_r) ...
0102             - apertureInfo.beam(i).posOfCornerBixel(1))/apertureInfo.bixelWidth;
0103         
0104         % find the bixel index that the leaves currently touch
0105         bixelIndLeftLeaf  = apertureInfo.beam(i).bixelIndMap((xPosIndLeftLeaf-1)*n+[1:n]');
0106         bixelIndRightLeaf = apertureInfo.beam(i).bixelIndMap((xPosIndRightLeaf-1)*n+[1:n]');
0107         
0108         if any(isnan(bixelIndLeftLeaf)) || any(isnan(bixelIndRightLeaf))
0109             error('cannot map leaf position to bixel index');
0110         end
0111         
0112         % store information in index vector for gradient calculation
0113         indVect(apertureInfo.beam(i).shape(j).vectorOffset+[1:n]-1) = bixelIndLeftLeaf;
0114         indVect(apertureInfo.beam(i).shape(j).vectorOffset+[1:n]-1+apertureInfo.totalNumOfLeafPairs) = bixelIndRightLeaf;
0115 
0116         % calculate opening fraction for every bixel in shape to construct
0117         % bixel weight vector
0118         
0119         coveredByLeftLeaf  = bsxfun(@minus,leftLeafPos,edges_l)  / updatedInfo.bixelWidth;
0120         coveredByRightLeaf = bsxfun(@minus,edges_r,rightLeafPos) / updatedInfo.bixelWidth;
0121 
0122         tempMap = 1 - (coveredByLeftLeaf  + abs(coveredByLeftLeaf))  / 2 ...
0123                     - (coveredByRightLeaf + abs(coveredByRightLeaf)) / 2;
0124             
0125         % find open bixels
0126         tempMapIx = tempMap > 0;
0127         
0128         currBixelIx = apertureInfo.beam(i).bixelIndMap(tempMapIx);
0129         w(currBixelIx) = w(currBixelIx) + tempMap(tempMapIx)*updatedInfo.beam(i).shape(j).weight;
0130     
0131         % save the tempMap (we need to apply a positivity operator !)
0132         updatedInfo.beam(i).shape(j).shapeMap = (tempMap  + abs(tempMap))  / 2;
0133         
0134         % increment shape index
0135         shapeInd = shapeInd +1;
0136     end
0137     
0138 end
0139 
0140 updatedInfo.bixelWeights = w;
0141 updatedInfo.bixelIndices = indVect;
0142 
0143 end
0144

| Generated by m2html © 2005