Skip to content

Commit

Permalink
fix bug in constructing GPbladeconstitutive
Browse files Browse the repository at this point in the history
  • Loading branch information
sean-engelstad committed Oct 8, 2024
1 parent d71c8ee commit e3be2cd
Show file tree
Hide file tree
Showing 5 changed files with 8 additions and 8 deletions.
4 changes: 2 additions & 2 deletions src/constitutive/TACSBladeStiffenedShellConstitutive.h
Original file line number Diff line number Diff line change
Expand Up @@ -1397,8 +1397,8 @@ class TACSBladeStiffenedShellConstitutive : public TACSShellConstitutive {
std::function<TacsScalar(const TacsScalar[],TacsScalar[])> evalLocalPanelBucklingStrainSens;
std::function<TacsScalar(const TacsScalar[],TacsScalar[])> evalGlobalPanelBucklingStrainSens;
std::function<TacsScalar(const TacsScalar[],TacsScalar[])> evalStiffenerCripplingStrainSens;
std::function<void(int,TacsScalar,const double[],const TacsScalar[],const TacsScalar[],int dvLen,TacsScalar[])> addLocalPanelBucklingDVSens;
std::function<void(int,TacsScalar,const double[],const TacsScalar[],const TacsScalar[],int dvLen,TacsScalar[])> addGlobalPanelBucklingDVSens;
std::function<void(int,TacsScalar,const double[],const TacsScalar[],const TacsScalar[],int,TacsScalar[])> addLocalPanelBucklingDVSens;
std::function<void(int,TacsScalar,const double[],const TacsScalar[],const TacsScalar[],int,TacsScalar[])> addGlobalPanelBucklingDVSens;
std::function<void(const TacsScalar, const TacsScalar[], TacsScalar[])> addStiffenerCripplingDVSens;

static const char* const constName; ///< Constitutive model name
Expand Down
7 changes: 2 additions & 5 deletions src/constitutive/TACSGPBladeStiffenedShellConstitutive.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -210,10 +210,7 @@ TacsScalar TACSGPBladeStiffenedShellConstitutive::_evalGlobalPanelBuckling(const
TacsScalar N12CritGlobal = computeCriticalGlobalShearLoad(
D11Global, D22p, b, rho0Global,
xiGlobal, gamma, zetaPanel
);

printf("N1CritGLobal = %.8e\n", N1CritGlobal);
printf("N12CritGlobal = %.8e\n", N12CritGlobal);
);

// compute the combined loading buckling failure index
return this->bucklingEnvelope(-panelStress[0], N1CritGlobal,
Expand Down Expand Up @@ -1439,7 +1436,7 @@ TACSGPBladeStiffenedShellConstitutive::computeStiffenerStiffnessRatioSens(
// backpropagate to intermediate quantities As, Is, zs, zn
TacsScalar AsSens = IstotSens * (zn - zs) * (zn - zs);
TacsScalar IsSens = IstotSens;
TacsScalar zsSens = IstotSens * As * 2.0 * (zs - zn);
TacsScalar zsSens = -1.0 * IstotSens * As * 2.0 * (zs - zn);
*znSens += IstotSens * As * 2.0 * (zn - zs);

// backpropgate from intermediate quantities to DVs or material, etc.
Expand Down
2 changes: 2 additions & 0 deletions tacs/constitutive.pyx
Original file line number Diff line number Diff line change
Expand Up @@ -1433,6 +1433,7 @@ cdef class GPBladeStiffenedShellConstitutive(ShellConstitutive):
int stiffenerThickNum = -1,
np.ndarray[int, ndim=1, mode='c'] stiffenerPlyFracNums = None,
int panelWidthNum = -1,
bool CPTstiffenerCrippling = False,
PanelGPs panelGPs = None,
):

Expand Down Expand Up @@ -1489,6 +1490,7 @@ cdef class GPBladeStiffenedShellConstitutive(ShellConstitutive):
panelWidth,
panelWidthNum,
flangeFraction,
CPTstiffenerCrippling,
panel_gp_ptr,
)
# copy pointers to all superclasses
Expand Down
1 change: 1 addition & 0 deletions tacs/cpp_headers/constitutive.pxd
Original file line number Diff line number Diff line change
Expand Up @@ -252,6 +252,7 @@ cdef extern from "TACSGPBladeStiffenedShellConstitutive.h":
TacsScalar, # panelWidth
int, # panelWidthNum
TacsScalar, # flangeFraction,
bool, # CPTstiffenerCrippling
TACSPanelGPs*, # panelGPs object container
)
TacsScalar nondimCriticalGlobalAxialLoad(TacsScalar rho_0, TacsScalar xi, TacsScalar gamma, TacsScalar zeta)
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -205,7 +205,7 @@ def get_con(self, ply):
self.stiffenerThicknessNum,
self.stiffenerPlyFracNums,
self.panelWidthNum,
self.panelGPs,
panelGPs=self.panelGPs,
)
# Set the KS weight really low so that all failure modes make a
# significant contribution to the failure function derivatives
Expand Down

0 comments on commit e3be2cd

Please sign in to comment.